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

#include <DisContField2D.h>

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

Public Member Functions

 DisContField2D ()
 DisContField2D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &graph2D, const std::string &variable, const bool SetUpJustDG=true, const bool DeclareCoeffPhysArrays=true)
 DisContField2D (const DisContField2D &In, const SpatialDomains::MeshGraphSharedPtr &graph2D, const std::string &variable, const bool SetUpJustDG=false, const bool DeclareCoeffPhysArrays=true)
 DisContField2D (const DisContField2D &In, const bool DeclareCoeffPhysArrays=true)
virtual ~DisContField2D ()
 Default destructor.
GlobalLinSysSharedPtr GetGlobalBndLinSys (const GlobalLinSysKey &mkey)
NekDouble L2_DGDeriv (const int dir, const Array< OneD, const NekDouble > &soln)
 Calculate the $ L^2 $ error of the $ Q_{\rm dir} $ derivative using the consistent DG evaluation of $ Q_{\rm dir} $.
void EvaluateHDGPostProcessing (Array< OneD, NekDouble > &outarray)
 Evaluate HDG post-processing to increase polynomial order of solution.
virtual ExpListSharedPtrv_GetTrace ()
- Public Member Functions inherited from Nektar::MultiRegions::ExpList2D
 ExpList2D ()
 Default constructor.
 ExpList2D (const ExpList2D &In, const bool DeclareCoeffPhysArrays=true)
 Copy constructor.
 ExpList2D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &graph2D, const bool DelcareCoeffPhysArrays=true, const std::string &var="DefaultVar")
 Sets up a list of local expansions based on an input mesh.
 ExpList2D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::ExpansionMap &expansions, const bool DeclareCoeffPhysArrays=true)
 Sets up a list of local expansions based on an expansion Map.
 ExpList2D (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &TriBa, const LibUtilities::BasisKey &TriBb, const LibUtilities::BasisKey &QuadBa, const LibUtilities::BasisKey &QuadBb, const SpatialDomains::MeshGraphSharedPtr &graph2D, const LibUtilities::PointsType TriNb=LibUtilities::SIZE_PointsType)
 Sets up a list of local expansions based on an input mesh and separately defined basiskeys.
 ExpList2D (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, const ExpListSharedPtr > &bndConstraint, const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > &bndCond, const LocalRegions::ExpansionVector &locexp, const SpatialDomains::MeshGraphSharedPtr &graph3D, const PeriodicMap &periodicFaces, const bool DeclareCoeffPhysArrays=true, const std::string variable="DefaultVar")
 ExpList2D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::CompositeMap &domain, const SpatialDomains::MeshGraphSharedPtr &graph3D, const std::string variable="DefaultVar")
 Specialised constructor for Neumann boundary conditions in DisContField3D and ContField3D.
virtual ~ExpList2D ()
 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::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 (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)
const Array< OneD, const
boost::shared_ptr< ExpList > > & 
GetBndCondExpansions ()
boost::shared_ptr< ExpList > & UpdateBndCondExpansion (int i)
void Upwind (const Array< OneD, const Array< OneD, NekDouble > > &Vec, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
void Upwind (const Array< OneD, const NekDouble > &Vn, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
boost::shared_ptr< ExpList > & GetTrace ()
boost::shared_ptr
< AssemblyMapDG > & 
GetTraceMap (void)
const Array< OneD, const int > & GetTraceBndMap (void)
void GetNormals (Array< OneD, Array< OneD, NekDouble > > &normals)
void AddTraceIntegral (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
void AddTraceIntegral (const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
void AddFwdBwdTraceIntegral (const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
void GetFwdBwdTracePhys (Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
void GetFwdBwdTracePhys (const Array< OneD, const NekDouble > &field, Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
void ExtractTracePhys (Array< OneD, NekDouble > &outarray)
void ExtractTracePhys (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
const Array< OneD, const
SpatialDomains::BoundaryConditionShPtr > & 
GetBndConditions ()
Array< OneD,
SpatialDomains::BoundaryConditionShPtr > & 
UpdateBndConditions ()
void EvaluateBoundaryConditions (const NekDouble time=0.0, const std::string varName="", const NekDouble=NekConstants::kNekUnsetDouble, const NekDouble=NekConstants::kNekUnsetDouble)
void GeneralMatrixOp (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
 This function calculates the result of the multiplication of a matrix of type specified by mkey with a vector given by inarray.
void GeneralMatrixOp_IterPerExp (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void SetUpPhysNormals ()
void GetBoundaryToElmtMap (Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
void GeneralGetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef, int NumHomoDir=0, Array< OneD, LibUtilities::BasisSharedPtr > &HomoBasis=LibUtilities::NullBasisSharedPtr1DArray, std::vector< NekDouble > &HomoLen=LibUtilities::NullNekDoubleVector, std::vector< unsigned int > &HomoZIDs=LibUtilities::NullUnsignedIntVector, std::vector< unsigned int > &HomoYIDs=LibUtilities::NullUnsignedIntVector)
const
NekOptimize::GlobalOptParamSharedPtr
GetGlobalOptParam (void)
map< int, RobinBCInfoSharedPtrGetRobinBCInfo ()
void GetPeriodicEntities (PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces=NullPeriodicMap)
std::vector
< LibUtilities::FieldDefinitionsSharedPtr
GetFieldDefinitions ()
void GetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
void AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
 Append the element data listed in elements fielddef->m_ElementIDs onto fielddata.
void AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, Array< OneD, NekDouble > &coeffs)
 Append the data in coeffs listed in elements fielddef->m_ElementIDs onto fielddata.
void ExtractElmtDataToCoeffs (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
 Extract the data in fielddata into the coeffs using the basic ExpList Elemental expansions rather than planes in homogeneous case.
void ExtractCoeffsToCoeffs (const boost::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
 Extract the data from fromField using fromExpList the coeffs using the basic ExpList Elemental expansions rather than planes in homogeneous case.
void ExtractDataToCoeffs (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
 Extract the data in fielddata into the coeffs.
boost::shared_ptr< ExpListGetSharedThisPtr ()
 Returns a shared pointer to the current object.
boost::shared_ptr
< LibUtilities::SessionReader
GetSession ()
 Returns the session object.
boost::shared_ptr
< LibUtilities::Comm
GetComm ()
 Returns the comm object.
SpatialDomains::MeshGraphSharedPtr GetGraph ()
LibUtilities::BasisSharedPtr GetHomogeneousBasis (void)
boost::shared_ptr< ExpList > & GetPlane (int n)

Protected Member Functions

void SetUpDG (const std::string="DefaultVar")
 Set up all DG member variables and maps.
bool SameTypeOfBoundaryConditions (const DisContField2D &In)
void GenerateBoundaryConditionExpansion (const SpatialDomains::MeshGraphSharedPtr &graph2D, const SpatialDomains::BoundaryConditions &bcs, const std::string &variable, const bool DeclareCoeffPhysArrays=true)
 This function discretises the boundary conditions by setting up a list of one-dimensional boundary expansions.
void FindPeriodicEdges (const SpatialDomains::BoundaryConditions &bcs, const std::string &variable)
 Determine the periodic edges and vertices for the given graph.
bool IsLeftAdjacentEdge (const int n, const int e)
virtual void v_GetFwdBwdTracePhys (const Array< OneD, const NekDouble > &field, 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 (Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
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)
 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_ExtractTracePhys (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This method extracts the trace (edges in 2D) from the field inarray and puts the values in outarray.
virtual void v_ExtractTracePhys (Array< OneD, NekDouble > &outarray)
virtual void v_FillBndCondFromField ()
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 > &EdgeID)
 Set up a list of element IDs and edge IDs that link to the boundary conditions.
virtual void v_GetPeriodicEntities (PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces)
 Obtain a copy of the periodic edges and vertices for this field.
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)
virtual map< int,
RobinBCInfoSharedPtr
v_GetRobinBCInfo ()
 Search through the edge expansions and identify which ones have Robin/Mixed type boundary conditions.
- Protected Member Functions inherited from Nektar::MultiRegions::ExpList2D
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)
 Upwind the Fwd and Bwd states based on the one- dimensional normal velocity field given by Vn.
void v_GetNormals (Array< OneD, Array< OneD, NekDouble > > &normals)
- 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 const Array< OneD,
const int > & 
v_GetTraceBndMap ()
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_LocalToGlobal (void)
virtual void v_GlobalToLocal (void)
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_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_SmoothField (Array< OneD, NekDouble > &field)
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 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_SetUpPhysNormals ()
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::ofstream &outfile, std::string var="")
virtual void v_WriteTecplotZone (std::ofstream &outfile, int expansion)
virtual void v_WriteTecplotField (std::ofstream &outfile, int expansion)
virtual void v_WriteTecplotConnectivity (std::ofstream &outfile, int expansion)
virtual void v_WriteVtkPieceHeader (std::ofstream &outfile, int expansion)
virtual void v_WriteVtkPieceData (std::ofstream &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

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
Array< OneD, Array< OneD,
unsigned int > > 
m_mapEdgeToElmn
Array< OneD, Array< OneD,
unsigned int > > 
m_signEdgeToElmn
Array< OneD,
StdRegions::Orientation
m_edgedir
std::set< int > m_boundaryEdges
 A set storing the global IDs of any boundary edges.
PeriodicMap m_periodicVerts
 A map which identifies groups of periodic vertices.
PeriodicMap m_periodicEdges
 A map which identifies pairs of periodic edges.
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
vector< bool > m_leftAdjacentEdges

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 DisContField2D.h.

Constructor & Destructor Documentation

Nektar::MultiRegions::DisContField2D::DisContField2D ( void  )
Nektar::MultiRegions::DisContField2D::DisContField2D ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr graph2D,
const std::string &  variable,
const bool  SetUpJustDG = true,
const bool  DeclareCoeffPhysArrays = true 
)

Definition at line 83 of file DisContField2D.cpp.

References Nektar::MultiRegions::ExpList::EvaluateBoundaryConditions(), FindPeriodicEdges(), GenerateBoundaryConditionExpansion(), Nektar::MultiRegions::ExpList::GetBoundaryToElmtMap(), m_bndCondExpansions, Nektar::MultiRegions::ExpList::m_session, SetUpDG(), and Nektar::MultiRegions::ExpList::SetUpPhysNormals().

: ExpList2D(pSession,graph2D,DeclareCoeffPhysArrays,variable),
{
if(variable.compare("DefaultVar") != 0) // do not set up BCs if default variable
{
GenerateBoundaryConditionExpansion(graph2D,bcs,variable,
DeclareCoeffPhysArrays);
if (DeclareCoeffPhysArrays)
{
}
// Find periodic edges for this variable.
FindPeriodicEdges(bcs, variable);
}
if (SetUpJustDG)
{
SetUpDG(variable);
}
else
{
// Set element edges to point to Robin BC edges if required.
int i, cnt;
Array<OneD, int> ElmtID, EdgeID;
GetBoundaryToElmtMap(ElmtID, EdgeID);
for(cnt = i = 0; i < m_bndCondExpansions.num_elements(); ++i)
{
int e;
locExpList = m_bndCondExpansions[i];
for(e = 0; e < locExpList->GetExpSize(); ++e)
{
(*m_exp)[ElmtID[cnt+e]]->
as<LocalRegions::Expansion2D>();
locExpList->GetExp(e)->
as<LocalRegions::Expansion1D>();
locExpList->GetExp(e)->
as<LocalRegions::Expansion> ();
exp2d->SetEdgeExp(EdgeID[cnt+e], exp);
exp1d->SetAdjacentElementExp(EdgeID[cnt+e], exp2d);
}
cnt += m_bndCondExpansions[i]->GetExpSize();
}
if(m_session->DefinesSolverInfo("PROJECTION"))
{
std::string ProjectStr =
m_session->GetSolverInfo("PROJECTION");
if((ProjectStr == "MixedCGDG") ||
(ProjectStr == "Mixed_CG_Discontinuous"))
{
}
else
{
}
}
else
{
}
}
}
Nektar::MultiRegions::DisContField2D::DisContField2D ( const DisContField2D In,
const SpatialDomains::MeshGraphSharedPtr graph2D,
const std::string &  variable,
const bool  SetUpJustDG = false,
const bool  DeclareCoeffPhysArrays = true 
)

Definition at line 173 of file DisContField2D.cpp.

References Nektar::MultiRegions::ExpList::EvaluateBoundaryConditions(), FindPeriodicEdges(), GenerateBoundaryConditionExpansion(), Nektar::MultiRegions::ExpList::GetBoundaryToElmtMap(), m_bndCondExpansions, m_boundaryEdges, m_globalBndMat, m_leftAdjacentEdges, m_periodicBwdCopy, m_periodicEdges, m_periodicFwdCopy, m_periodicVerts, Nektar::MultiRegions::ExpList::m_session, m_trace, m_traceMap, SameTypeOfBoundaryConditions(), SetUpDG(), and Nektar::MultiRegions::ExpList::SetUpPhysNormals().

: ExpList2D(In,DeclareCoeffPhysArrays),
{
// Set up boundary conditions for this variable.
// Do not set up BCs if default variable
if(variable.compare("DefaultVar") != 0)
{
GenerateBoundaryConditionExpansion(graph2D, bcs, variable);
if (DeclareCoeffPhysArrays)
{
}
{
// Find periodic edges for this variable.
FindPeriodicEdges(bcs, variable);
if(SetUpJustDG)
{
}
else
{
// set elmt edges to point to robin bc edges if required.
int i, cnt = 0;
Array<OneD, int> ElmtID,EdgeID;
GetBoundaryToElmtMap(ElmtID,EdgeID);
for(i = 0; i < m_bndCondExpansions.num_elements(); ++i)
{
int e;
locExpList = m_bndCondExpansions[i];
for(e = 0; e < locExpList->GetExpSize(); ++e)
{
= (*m_exp)[ElmtID[cnt+e]]->
as<LocalRegions::Expansion2D>();
= locExpList->GetExp(e)->
as<LocalRegions::Expansion1D>();
= locExpList->GetExp(e)->
as<LocalRegions::Expansion> ();
exp2d->SetEdgeExp(EdgeID[cnt+e],exp);
exp1d->SetAdjacentElementExp(EdgeID[cnt+e],exp2d);
}
cnt += m_bndCondExpansions[i]->GetExpSize();
}
if(m_session->DefinesSolverInfo("PROJECTION"))
{
std::string ProjectStr =
m_session->GetSolverInfo("PROJECTION");
if((ProjectStr == "MixedCGDG") ||
(ProjectStr == "Mixed_CG_Discontinuous"))
{
}
else
{
}
}
else
{
}
}
}
else
{
if(SetUpJustDG)
{
m_globalBndMat = In.m_globalBndMat;
m_trace = In.m_trace;
m_traceMap = In.m_traceMap;
m_periodicEdges = In.m_periodicEdges;
m_periodicVerts = In.m_periodicVerts;
m_periodicFwdCopy = In.m_periodicFwdCopy;
m_periodicBwdCopy = In.m_periodicBwdCopy;
m_boundaryEdges = In.m_boundaryEdges;
m_leftAdjacentEdges = In.m_leftAdjacentEdges;
}
else
{
m_globalBndMat = In.m_globalBndMat;
m_trace = In.m_trace;
m_traceMap = In.m_traceMap;
m_periodicEdges = In.m_periodicEdges;
m_periodicVerts = In.m_periodicVerts;
m_periodicFwdCopy = In.m_periodicFwdCopy;
m_periodicBwdCopy = In.m_periodicBwdCopy;
m_boundaryEdges = In.m_boundaryEdges;
m_leftAdjacentEdges = In.m_leftAdjacentEdges;
// set elmt edges to point to robin bc edges if required.
int i, cnt = 0;
Array<OneD, int> ElmtID, EdgeID;
GetBoundaryToElmtMap(ElmtID, EdgeID);
for(i = 0; i < m_bndCondExpansions.num_elements(); ++i)
{
int e;
locExpList = m_bndCondExpansions[i];
for(e = 0; e < locExpList->GetExpSize(); ++e)
{
= (*m_exp)[ElmtID[cnt+e]]->
as<LocalRegions::Expansion2D>();
= locExpList->GetExp(e)->
as<LocalRegions::Expansion1D>();
= locExpList->GetExp(e)->
as<LocalRegions::Expansion> ();
exp2d->SetEdgeExp(EdgeID[cnt+e],exp);
exp1d->SetAdjacentElementExp(EdgeID[cnt+e],exp2d);
}
cnt += m_bndCondExpansions[i]->GetExpSize();
}
}
}
}
}
Nektar::MultiRegions::DisContField2D::DisContField2D ( const DisContField2D In,
const bool  DeclareCoeffPhysArrays = true 
)

Definition at line 60 of file DisContField2D.cpp.

References m_trace.

: ExpList2D (In,DeclareCoeffPhysArrays),
m_bndCondExpansions (In.m_bndCondExpansions),
m_bndConditions (In.m_bndConditions),
m_globalBndMat (In.m_globalBndMat),
m_traceMap (In.m_traceMap),
m_boundaryEdges (In.m_boundaryEdges),
m_periodicVerts (In.m_periodicVerts),
m_periodicEdges (In.m_periodicEdges),
m_periodicFwdCopy (In.m_periodicFwdCopy),
m_periodicBwdCopy (In.m_periodicBwdCopy),
m_leftAdjacentEdges (In.m_leftAdjacentEdges)
{
if (In.m_trace)
{
*boost::dynamic_pointer_cast<ExpList1D>(In.m_trace),
DeclareCoeffPhysArrays);
}
}
Nektar::MultiRegions::DisContField2D::~DisContField2D ( )
virtual

Default destructor.

Definition at line 322 of file DisContField2D.cpp.

{
}

Member Function Documentation

void Nektar::MultiRegions::DisContField2D::EvaluateHDGPostProcessing ( Array< OneD, NekDouble > &  outarray)

Evaluate HDG post-processing to increase polynomial order of solution.

This function takes the solution (assumed to be one order lower) in physical space, and postprocesses at the current polynomial order by solving the system:

\[ \begin{aligned} (\nabla w, \nabla u^*) &= (\nabla w, u), \\ \langle \nabla u^*, 1 \rangle &= \langle \nabla u, 1 \rangle \end{aligned} \]

where $ u $ corresponds with the current solution as stored inside m_coeffs.

Parameters
outarrayThe resulting field $ u^* $.

Definition at line 1965 of file DisContField2D.cpp.

References ASSERTL0, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::StdRegions::eInvLaplacianWithUnityMean, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrtho_B, Nektar::LibUtilities::eQuadrilateral, Nektar::LibUtilities::eTriangle, Nektar::eWrapper, Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::NekVector< DataType >::GetPtr(), Nektar::MultiRegions::ExpList::m_coeff_offset, Nektar::MultiRegions::ExpList::m_coeffs, Nektar::MultiRegions::ExpList::m_offset_elmt_id, m_trace, m_traceMap, Vmath::Vadd(), and Vmath::Vcopy().

{
int i,cnt,e,ncoeff_edge;
Array<OneD, NekDouble> force, out_tmp, qrhs, qrhs1;
Array<OneD, Array< OneD, StdRegions::StdExpansionSharedPtr> >
&elmtToTrace = m_traceMap->GetElmtToTrace();
int eid,nq_elmt, nm_elmt;
int LocBndCoeffs = m_traceMap->GetNumLocalBndCoeffs();
Array<OneD, NekDouble> loc_lambda(LocBndCoeffs), edge_lambda;
Array<OneD, NekDouble> tmp_coeffs;
m_traceMap->GlobalToLocalBnd(m_trace->GetCoeffs(),loc_lambda);
edge_lambda = loc_lambda;
// Calculate Q using standard DG formulation.
for(i = cnt = 0; i < GetExpSize(); ++i)
{
eid = m_offset_elmt_id[i];
nq_elmt = (*m_exp)[eid]->GetTotPoints();
nm_elmt = (*m_exp)[eid]->GetNcoeffs();
qrhs = Array<OneD, NekDouble>(nq_elmt);
qrhs1 = Array<OneD, NekDouble>(nq_elmt);
force = Array<OneD, NekDouble>(2*nm_elmt);
out_tmp = force + nm_elmt;
int num_points0 = (*m_exp)[eid]->GetBasis(0)->GetNumPoints();
int num_points1 = (*m_exp)[eid]->GetBasis(1)->GetNumPoints();
int num_modes0 = (*m_exp)[eid]->GetBasis(0)->GetNumModes();
int num_modes1 = (*m_exp)[eid]->GetBasis(1)->GetNumModes();
// Probably a better way of setting up lambda than this. Note
// cannot use PutCoeffsInToElmts since lambda space is mapped
// during the solve.
int nEdges = (*m_exp)[i]->GetNedges();
Array<OneD, Array<OneD, NekDouble> > edgeCoeffs(nEdges);
for(e = 0; e < (*m_exp)[eid]->GetNedges(); ++e)
{
edgedir = (*m_exp)[eid]->GetEorient(e);
ncoeff_edge = elmtToTrace[eid][e]->GetNcoeffs();
edgeCoeffs[e] = Array<OneD, NekDouble>(ncoeff_edge);
Vmath::Vcopy(ncoeff_edge, edge_lambda, 1, edgeCoeffs[e], 1);
elmtToTrace[eid][e]->SetCoeffsToOrientation(
edgedir, edgeCoeffs[e], edgeCoeffs[e]);
edge_lambda = edge_lambda + ncoeff_edge;
}
//creating orthogonal expansion (checking if we have quads or triangles)
LibUtilities::ShapeType shape = (*m_exp)[eid]->DetShapeType();
switch(shape)
{
{
LibUtilities::BasisKey BkeyQ1(LibUtilities::eOrtho_A, num_modes0, PkeyQ1);
LibUtilities::BasisKey BkeyQ2(LibUtilities::eOrtho_A, num_modes1, PkeyQ2);
SpatialDomains::QuadGeomSharedPtr qGeom = boost::dynamic_pointer_cast<SpatialDomains::QuadGeom>((*m_exp)[eid]->GetGeom());
}
break;
{
LibUtilities::BasisKey BkeyT1(LibUtilities::eOrtho_A, num_modes0, PkeyT1);
LibUtilities::BasisKey BkeyT2(LibUtilities::eOrtho_B, num_modes1, PkeyT2);
SpatialDomains::TriGeomSharedPtr tGeom = boost::dynamic_pointer_cast<SpatialDomains::TriGeom>((*m_exp)[eid]->GetGeom());
}
break;
default:
ASSERTL0(false, "Wrong shape type, HDG postprocessing is not implemented");
};
//SpatialDomains::QuadGeomSharedPtr qGeom = boost::dynamic_pointer_cast<SpatialDomains::QuadGeom>((*m_exp)[eid]->GetGeom());
//LocalRegions::QuadExpSharedPtr ppExp =
// MemoryManager<LocalRegions::QuadExp>::AllocateSharedPtr(BkeyQ1, BkeyQ2, qGeom);
//Orthogonal expansion created
//In case lambdas are causing the trouble, try PhysDeriv instead of DGDeriv
//===============================================================================================
//(*m_exp)[eid]->BwdTrans(tmp_coeffs = m_coeffs + m_coeff_offset[eid],(*m_exp)[eid]->UpdatePhys());
//(*m_exp)[eid]->PhysDeriv((*m_exp)[eid]->GetPhys(), qrhs, qrhs1);
//ppExp->IProductWRTDerivBase(0,qrhs,force);
//ppExp->IProductWRTDerivBase(1,qrhs1,out_tmp);
//===============================================================================================
//DGDeriv
// (d/dx w, d/dx q_0)
(*m_exp)[eid]->DGDeriv(
0,tmp_coeffs = m_coeffs + m_coeff_offset[eid],
elmtToTrace[eid], edgeCoeffs, out_tmp);
(*m_exp)[eid]->BwdTrans(out_tmp,qrhs);
//(*m_exp)[eid]->IProductWRTDerivBase(0,qrhs,force);
ppExp->IProductWRTDerivBase(0,qrhs,force);
// + (d/dy w, d/dy q_1)
(*m_exp)[eid]->DGDeriv(
1,tmp_coeffs = m_coeffs + m_coeff_offset[eid],
elmtToTrace[eid], edgeCoeffs, out_tmp);
(*m_exp)[eid]->BwdTrans(out_tmp,qrhs);
//(*m_exp)[eid]->IProductWRTDerivBase(1,qrhs,out_tmp);
ppExp->IProductWRTDerivBase(1,qrhs,out_tmp);
Vmath::Vadd(nm_elmt,force,1,out_tmp,1,force,1);
// determine force[0] = (1,u)
(*m_exp)[eid]->BwdTrans(
tmp_coeffs = m_coeffs + m_coeff_offset[eid],qrhs);
force[0] = (*m_exp)[eid]->Integral(qrhs);
// multiply by inverse Laplacian matrix
// get matrix inverse
LocalRegions::MatrixKey lapkey(StdRegions::eInvLaplacianWithUnityMean, ppExp->DetShapeType(), *ppExp);
DNekScalMatSharedPtr lapsys = ppExp->GetLocMatrix(lapkey);
NekVector<NekDouble> in (nm_elmt,force,eWrapper);
NekVector<NekDouble> out(nm_elmt);
out = (*lapsys)*in;
// Transforming back to modified basis
Array<OneD, NekDouble> work(nq_elmt);
ppExp->BwdTrans(out.GetPtr(), work);
(*m_exp)[eid]->FwdTrans(work, tmp_coeffs = outarray + m_coeff_offset[eid]);
}
}
void Nektar::MultiRegions::DisContField2D::FindPeriodicEdges ( const SpatialDomains::BoundaryConditions bcs,
const std::string &  variable 
)
protected

Determine the periodic edges and vertices for the given graph.

Note that much of this routine is the same as the three-dimensional version, which therefore has much better documentation.

Parameters
bcsInformation about the boundary conditions.
variableSpecifies the field.
See Also
DisContField3D::FindPeriodicFaces

Definition at line 679 of file DisContField2D.cpp.

References ASSERTL0, Nektar::StdRegions::eBackwards, Nektar::StdRegions::eForwards, Nektar::StdRegions::eNoOrientation, Nektar::SpatialDomains::ePeriodic, Nektar::MultiRegions::ExpList::GetBoundaryCondition(), Nektar::SpatialDomains::BoundaryConditions::GetBoundaryConditions(), Nektar::SpatialDomains::BoundaryConditions::GetBoundaryRegions(), Nektar::iterator, Nektar::MultiRegions::ExpList::m_graph, m_periodicEdges, m_periodicVerts, Nektar::MultiRegions::ExpList::m_session, Nektar::LibUtilities::ReduceSum, and Vmath::Vsum().

Referenced by DisContField2D().

{
= boost::dynamic_pointer_cast<
SpatialDomains::BoundaryRegionCollection::const_iterator it;
m_session->GetComm()->GetRowComm();
m_session->GetCompositeOrdering();
m_session->GetBndRegionOrdering();
m_graph->GetComposites();
// Unique collection of pairs of periodic composites (i.e. if
// composites 1 and 2 are periodic then this map will contain either
// the pair (1,2) or (2,1) but not both).
map<int,int> perComps;
map<int,vector<int> > allVerts;
set<int> locVerts;
map<int,StdRegions::Orientation> allEdges;
int region1ID, region2ID, i, j, k, cnt;
// Set up a set of all local verts and edges.
for(i = 0; i < (*m_exp).size(); ++i)
{
for(j = 0; j < (*m_exp)[i]->GetNverts(); ++j)
{
int id = (*m_exp)[i]->GetGeom()->GetVid(j);
locVerts.insert(id);
}
}
// Construct list of all periodic pairs local to this process.
for (it = bregions.begin(); it != bregions.end(); ++it)
{
bconditions, it->first, variable);
if (locBCond->GetBoundaryConditionType()
{
continue;
}
// Identify periodic boundary region IDs.
region1ID = it->first;
region2ID = boost::static_pointer_cast<
locBCond)->m_connectedBoundaryRegion;
// From this identify composites. Note that in serial this will
// be an empty map.
int cId1, cId2;
if (vComm->GetSize() == 1)
{
cId1 = it->second->begin()->first;
cId2 = bregions.find(region2ID)->second->begin()->first;
}
else
{
cId1 = bndRegOrder.find(region1ID)->second[0];
cId2 = bndRegOrder.find(region2ID)->second[0];
}
ASSERTL0(it->second->size() == 1,
"Boundary region "+boost::lexical_cast<string>(
region1ID)+" should only contain 1 composite.");
// Construct set containing all periodic edges on this process
SpatialDomains::Composite c = it->second->begin()->second;
vector<unsigned int> tmpOrder;
for (i = 0; i < c->size(); ++i)
{
boost::dynamic_pointer_cast<
ASSERTL0(segGeom, "Unable to cast to shared ptr");
graph2D->GetElementsFromEdge(segGeom);
ASSERTL0(elmt->size() == 1,
"The periodic boundaries belong to "
"more than one element of the mesh");
boost::dynamic_pointer_cast<SpatialDomains::Geometry2D>(
(*elmt)[0]->m_Element);
allEdges[(*c)[i]->GetGlobalID()] =
geom->GetEorient((*elmt)[0]->m_EdgeIndx);
// In serial mesh partitioning will not have occurred so
// need to fill composite ordering map manually.
if (vComm->GetSize() == 1)
{
tmpOrder.push_back((*c)[i]->GetGlobalID());
}
vector<int> vertList(2);
vertList[0] = segGeom->GetVid(0);
vertList[1] = segGeom->GetVid(1);
allVerts[(*c)[i]->GetGlobalID()] = vertList;
}
if (vComm->GetSize() == 1)
{
compOrder[it->second->begin()->first] = tmpOrder;
}
// See if we already have either region1 or region2 stored in
// perComps map.
if (perComps.count(cId1) == 0)
{
if (perComps.count(cId2) == 0)
{
perComps[cId1] = cId2;
}
else
{
std::stringstream ss;
ss << "Boundary region " << cId2 << " should be "
<< "periodic with " << perComps[cId2] << " but "
<< "found " << cId1 << " instead!";
ASSERTL0(perComps[cId2] == cId1, ss.str());
}
}
else
{
std::stringstream ss;
ss << "Boundary region " << cId1 << " should be "
<< "periodic with " << perComps[cId1] << " but "
<< "found " << cId2 << " instead!";
ASSERTL0(perComps[cId1] == cId1, ss.str());
}
}
// Process local edge list to obtain relative edge orientations.
int n = vComm->GetSize();
int p = vComm->GetRank();
int totEdges;
Array<OneD, int> edgecounts(n,0);
Array<OneD, int> edgeoffset(n,0);
Array<OneD, int> vertoffset(n,0);
edgecounts[p] = allEdges.size();
vComm->AllReduce(edgecounts, LibUtilities::ReduceSum);
edgeoffset[0] = 0;
for (i = 1; i < n; ++i)
{
edgeoffset[i] = edgeoffset[i-1] + edgecounts[i-1];
}
totEdges = Vmath::Vsum(n, edgecounts, 1);
Array<OneD, int> edgeIds (totEdges, 0);
Array<OneD, int> edgeOrient(totEdges, 0);
Array<OneD, int> edgeVerts (totEdges, 0);
for (i = 0, sIt = allEdges.begin(); sIt != allEdges.end(); ++sIt)
{
edgeIds [edgeoffset[p] + i ] = sIt->first;
edgeOrient[edgeoffset[p] + i ] = sIt->second;
edgeVerts [edgeoffset[p] + i++] = allVerts[sIt->first].size();
}
vComm->AllReduce(edgeIds, LibUtilities::ReduceSum);
vComm->AllReduce(edgeOrient, LibUtilities::ReduceSum);
vComm->AllReduce(edgeVerts, LibUtilities::ReduceSum);
// Calculate number of vertices on each processor.
Array<OneD, int> procVerts(n,0);
int nTotVerts;
// Note if there are no periodic edges at all calling Vsum will
// cause a segfault.
if (totEdges > 0)
{
nTotVerts = Vmath::Vsum(totEdges, edgeVerts, 1);
}
else
{
nTotVerts = 0;
}
for (i = 0; i < n; ++i)
{
if (edgecounts[i] > 0)
{
procVerts[i] = Vmath::Vsum(
edgecounts[i], edgeVerts + edgeoffset[i], 1);
}
else
{
procVerts[i] = 0;
}
}
vertoffset[0] = 0;
for (i = 1; i < n; ++i)
{
vertoffset[i] = vertoffset[i-1] + procVerts[i-1];
}
Array<OneD, int> vertIds(nTotVerts, 0);
for (i = 0, sIt = allEdges.begin(); sIt != allEdges.end(); ++sIt)
{
for (j = 0; j < allVerts[sIt->first].size(); ++j)
{
vertIds[vertoffset[p] + i++] = allVerts[sIt->first][j];
}
}
vComm->AllReduce(vertIds, LibUtilities::ReduceSum);
// For simplicity's sake create a map of edge id -> orientation.
map<int, StdRegions::Orientation> orientMap;
map<int, vector<int> > vertMap;
for (cnt = i = 0; i < totEdges; ++i)
{
ASSERTL0(orientMap.count(edgeIds[i]) == 0,
"Already found edge in orientation map!");
orientMap[edgeIds[i]] = (StdRegions::Orientation)edgeOrient[i];
vector<int> verts(edgeVerts[i]);
for (j = 0; j < edgeVerts[i]; ++j)
{
verts[j] = vertIds[cnt++];
}
vertMap[edgeIds[i]] = verts;
}
// Go through list of composites and figure out which edges are
// parallel from original ordering in session file. This includes
// composites which are not necessarily on this process.
map<int,int>::const_iterator oIt;
map<int,int> allCompPairs;
// Store temporary map of periodic vertices which will hold all
// periodic vertices on the entire mesh so that doubly periodic
// vertices can be counted properly across partitions. Local
// vertices are copied into m_periodicVerts at the end of the
// function.
PeriodicMap periodicVerts;
for (cIt = perComps.begin(); cIt != perComps.end(); ++cIt)
{
const int id1 = cIt->first;
const int id2 = cIt->second;
std::string id1s = boost::lexical_cast<string>(id1);
std::string id2s = boost::lexical_cast<string>(id2);
if (compMap.count(id1) > 0)
{
c[0] = compMap[id1];
}
if (compMap.count(id2) > 0)
{
c[1] = compMap[id2];
}
ASSERTL0(c[0] || c[1],
"Both composites not found on this process!");
// Loop over composite ordering to construct list of all
// periodic edges regardless of whether they are on this
// process.
map<int,int> compPairs;
ASSERTL0(compOrder.count(id1) > 0,
"Unable to find composite "+id1s+" in order map.");
ASSERTL0(compOrder.count(id2) > 0,
"Unable to find composite "+id2s+" in order map.");
ASSERTL0(compOrder[id1].size() == compOrder[id2].size(),
"Periodic composites "+id1s+" and "+id2s+
" should have the same number of elements.");
ASSERTL0(compOrder[id1].size() > 0,
"Periodic composites "+id1s+" and "+id2s+
" are empty!");
// TODO: Add more checks.
for (i = 0; i < compOrder[id1].size(); ++i)
{
int eId1 = compOrder[id1][i];
int eId2 = compOrder[id2][i];
ASSERTL0(compPairs.count(eId1) == 0,
"Already paired.");
if (compPairs.count(eId2) != 0)
{
ASSERTL0(compPairs[eId2] == eId1, "Pairing incorrect");
}
compPairs[eId1] = eId2;
}
// Construct set of all edges that we have locally on this
// processor.
set<int> locEdges;
for (i = 0; i < 2; ++i)
{
if (!c[i])
{
continue;
}
if (c[i]->size() > 0)
{
for (j = 0; j < c[i]->size(); ++j)
{
locEdges.insert(c[i]->at(j)->GetGlobalID());
}
}
}
// Loop over all edges in the geometry composite.
for (pIt = compPairs.begin(); pIt != compPairs.end(); ++pIt)
{
int ids [2] = {pIt->first, pIt->second};
bool local[2] = {locEdges.count(pIt->first) > 0,
locEdges.count(pIt->second) > 0};
ASSERTL0(orientMap.count(ids[0]) > 0 &&
orientMap.count(ids[1]) > 0,
"Unable to find edge in orientation map.");
allCompPairs[pIt->first ] = pIt->second;
allCompPairs[pIt->second] = pIt->first;
for (i = 0; i < 2; ++i)
{
if (!local[i])
{
continue;
}
int other = (i+1) % 2;
orientMap[ids[i]] == orientMap[ids[other]] ?
PeriodicEntity ent(ids [other], o,
local[other]);
m_periodicEdges[ids[i]].push_back(ent);
}
for (i = 0; i < 2; ++i)
{
int other = (i+1) % 2;
orientMap[ids[i]] == orientMap[ids[other]] ?
// Determine periodic vertices.
vector<int> perVerts1 = vertMap[ids[i]];
vector<int> perVerts2 = vertMap[ids[other]];
map<int, pair<int, bool> > tmpMap;
map<int, pair<int, bool> >::iterator mIt;
{
tmpMap[perVerts1[0]] = make_pair(
perVerts2[0], locVerts.count(perVerts2[0]) > 0);
tmpMap[perVerts1[1]] = make_pair(
perVerts2[1], locVerts.count(perVerts2[1]) > 0);
}
else
{
tmpMap[perVerts1[0]] = make_pair(
perVerts2[1], locVerts.count(perVerts2[1]) > 0);
tmpMap[perVerts1[1]] = make_pair(
perVerts2[0], locVerts.count(perVerts2[0]) > 0);
}
for (mIt = tmpMap.begin(); mIt != tmpMap.end(); ++mIt)
{
// See if this vertex has been recorded already.
PeriodicEntity ent2(mIt->second.first,
mIt->second.second);
PeriodicMap::iterator perIt = periodicVerts.find(
mIt->first);
if (perIt == periodicVerts.end())
{
periodicVerts[mIt->first].push_back(ent2);
perIt = periodicVerts.find(mIt->first);
}
else
{
bool doAdd = true;
for (j = 0; j < perIt->second.size(); ++j)
{
if (perIt->second[j].id == mIt->second.first)
{
doAdd = false;
break;
}
}
if (doAdd)
{
perIt->second.push_back(ent2);
}
}
}
}
}
}
Array<OneD, int> pairSizes(n, 0);
pairSizes[p] = allCompPairs.size();
vComm->AllReduce(pairSizes, LibUtilities::ReduceSum);
int totPairSizes = Vmath::Vsum(n, pairSizes, 1);
Array<OneD, int> pairOffsets(n, 0);
pairOffsets[0] = 0;
for (i = 1; i < n; ++i)
{
pairOffsets[i] = pairOffsets[i-1] + pairSizes[i-1];
}
Array<OneD, int> first (totPairSizes, 0);
Array<OneD, int> second(totPairSizes, 0);
cnt = pairOffsets[p];
for (pIt = allCompPairs.begin(); pIt != allCompPairs.end(); ++pIt)
{
first [cnt ] = pIt->first;
second[cnt++] = pIt->second;
}
vComm->AllReduce(first, LibUtilities::ReduceSum);
vComm->AllReduce(second, LibUtilities::ReduceSum);
allCompPairs.clear();
for(cnt = 0; cnt < totPairSizes; ++cnt)
{
allCompPairs[first[cnt]] = second[cnt];
}
// Search for periodic vertices and edges which are not in a
// periodic composite but lie in this process. First, loop over all
// information we have from other processors.
for (cnt = i = 0; i < totEdges; ++i)
{
int edgeId = edgeIds[i];
ASSERTL0(allCompPairs.count(edgeId) > 0,
"Unable to find matching periodic edge.");
int perEdgeId = allCompPairs[edgeId];
for (j = 0; j < edgeVerts[i]; ++j, ++cnt)
{
int vId = vertIds[cnt];
PeriodicMap::iterator perId = periodicVerts.find(vId);
if (perId == periodicVerts.end())
{
// This vertex is not included in the map. Figure out
// which vertex it is supposed to be periodic
// with. perEdgeId is the edge ID which is periodic with
// edgeId. The logic is much the same as the loop above.
int perVertexId =
orientMap[edgeId] == orientMap[perEdgeId] ?
vertMap[perEdgeId][(j+1)%2] : vertMap[perEdgeId][j];
PeriodicEntity ent(perVertexId,
locVerts.count(perVertexId) > 0);
periodicVerts[vId].push_back(ent);
}
}
}
// Loop over all periodic vertices to complete connectivity
// information.
PeriodicMap::iterator perIt, perIt2;
for (perIt = periodicVerts.begin();
perIt != periodicVerts.end(); ++perIt)
{
// Loop over associated vertices.
for (i = 0; i < perIt->second.size(); ++i)
{
perIt2 = periodicVerts.find(perIt->second[i].id);
ASSERTL0(perIt2 != periodicVerts.end(),
"Couldn't find periodic vertex.");
for (j = 0; j < perIt2->second.size(); ++j)
{
if (perIt2->second[j].id == perIt->first)
{
continue;
}
bool doAdd = true;
for (k = 0; k < perIt->second.size(); ++k)
{
if (perIt2->second[j].id == perIt->second[k].id)
{
doAdd = false;
break;
}
}
if (doAdd)
{
perIt->second.push_back(perIt2->second[j]);
}
}
}
}
// Do one final loop over periodic vertices to remove non-local
// vertices from map.
for (perIt = periodicVerts.begin();
perIt != periodicVerts.end(); ++perIt)
{
if (locVerts.count(perIt->first) > 0)
{
m_periodicVerts.insert(*perIt);
}
}
}
void Nektar::MultiRegions::DisContField2D::GenerateBoundaryConditionExpansion ( const SpatialDomains::MeshGraphSharedPtr graph2D,
const SpatialDomains::BoundaryConditions bcs,
const std::string &  variable,
const bool  DeclareCoeffPhysArrays = true 
)
protected

This function discretises the boundary conditions by setting up a list of one-dimensional boundary expansions.

According to their boundary region, the separate segmental boundary expansions are bundled together in an object of the class MultiRegions::ExpList1D.

Parameters
graph2DA mesh, containing information about the domain and the spectral/hp element expansion.
bcsAn entity containing information about the boundaries and boundary conditions.
variableAn optional parameter to indicate for which variable the boundary conditions should be discretised.

Definition at line 604 of file DisContField2D.cpp.

References Nektar::SpatialDomains::eCalcBC, Nektar::SpatialDomains::eDirichlet, Nektar::SpatialDomains::eI, Nektar::SpatialDomains::ePeriodic, Nektar::MultiRegions::ExpList::GetBoundaryCondition(), Nektar::SpatialDomains::BoundaryConditions::GetBoundaryConditions(), Nektar::SpatialDomains::BoundaryConditions::GetBoundaryRegions(), m_bndCondExpansions, m_bndConditions, and Nektar::MultiRegions::ExpList::SetUpPhysNormals().

Referenced by DisContField2D().

{
int cnt = 0;
SpatialDomains::BoundaryRegionCollection::const_iterator it;
// count the number of non-periodic boundary regions
for (it = bregions.begin(); it != bregions.end(); ++it)
{
bc = GetBoundaryCondition(bconditions, it->first, variable);
if (bc->GetBoundaryConditionType() != SpatialDomains::ePeriodic)
{
cnt++;
}
}
Array<OneD, MultiRegions::ExpListSharedPtr>(cnt);
Array<OneD, SpatialDomains::BoundaryConditionShPtr>(cnt);
cnt = 0;
// list non-periodic boundaries
for (it = bregions.begin(); it != bregions.end(); ++it)
{
bc = GetBoundaryCondition(bconditions, it->first, variable);
if (bc->GetBoundaryConditionType() != SpatialDomains::ePeriodic)
{
::AllocateSharedPtr(*(it->second), graph2D,
DeclareCoeffPhysArrays, variable);
// Set up normals on non-Dirichlet boundary conditions
if(bc->GetBoundaryConditionType() !=
{
}
m_bndCondExpansions[cnt] = locExpList;
m_bndConditions[cnt] = bc;
m_bndConditions[cnt++]->GetUserDefined();
if (type == SpatialDomains::eI ||
{
}
}
}
}
GlobalLinSysSharedPtr Nektar::MultiRegions::DisContField2D::GetGlobalBndLinSys ( const GlobalLinSysKey mkey)

Definition at line 326 of file DisContField2D.cpp.

References ASSERTL0, ASSERTL1, Nektar::StdRegions::eHybridDGHelmBndLam, Nektar::MultiRegions::ExpList::GenGlobalBndLinSys(), Nektar::MultiRegions::GlobalLinSysKey::GetGlobalSysSolnType(), Nektar::MultiRegions::GlobalMatrixKey::GetMatrixType(), Nektar::iterator, m_globalBndMat, and m_traceMap.

Referenced by v_HelmSolve().

{
ASSERTL0(mkey.GetMatrixType() == StdRegions::eHybridDGHelmBndLam,
"Routine currently only tested for HybridDGHelmholtz");
ASSERTL1(mkey.GetGlobalSysSolnType() ==
m_traceMap->GetGlobalSysSolnType(),
"The local to global map is not set up for the requested "
"solution type");
GlobalLinSysMap::iterator matrixIter = m_globalBndMat->find(mkey);
if(matrixIter == m_globalBndMat->end())
{
glo_matrix = GenGlobalBndLinSys(mkey,m_traceMap);
(*m_globalBndMat)[mkey] = glo_matrix;
}
else
{
glo_matrix = matrixIter->second;
}
return glo_matrix;
}
bool Nektar::MultiRegions::DisContField2D::IsLeftAdjacentEdge ( const int  n,
const int  e 
)
protected

Definition at line 1236 of file DisContField2D.cpp.

References ASSERTL2, Nektar::iterator, m_boundaryEdges, Nektar::MultiRegions::ExpList::m_exp, m_periodicEdges, m_trace, and m_traceMap.

Referenced by SetUpDG(), and v_AddFwdBwdTraceIntegral().

{
m_traceMap->GetElmtToTrace()[n][e]->
as<LocalRegions::Expansion1D>();
bool fwd = true;
if (traceEl->GetLeftAdjacentElementEdge () == -1 ||
traceEl->GetRightAdjacentElementEdge() == -1)
{
// Boundary edge (1 connected element). Do nothing in
// serial.
it = m_boundaryEdges.find(traceEl->GetElmtId());
// If the edge does not have a boundary condition set on
// it, then assume it is a partition edge.
if (it == m_boundaryEdges.end())
{
int traceGeomId = traceEl->GetGeom1D()->GetGlobalID();
traceGeomId);
if (pIt != m_periodicEdges.end() && !pIt->second[0].isLocal)
{
fwd = traceGeomId == min(traceGeomId,pIt->second[0].id);
}
else
{
int offset = m_trace->GetPhys_Offset(traceEl->GetElmtId());
fwd = m_traceMap->
GetTraceToUniversalMapUnique(offset) >= 0;
}
}
}
else if (traceEl->GetLeftAdjacentElementEdge () != -1 &&
traceEl->GetRightAdjacentElementEdge() != -1)
{
// Non-boundary edge (2 connected elements).
fwd = dynamic_cast<Nektar::StdRegions::StdExpansion*>
(traceEl->GetLeftAdjacentElementExp().get()) ==
(*m_exp)[n].get();
}
else
{
ASSERTL2(false, "Unconnected trace element!");
}
return fwd;
}
NekDouble Nektar::MultiRegions::DisContField2D::L2_DGDeriv ( const int  dir,
const Array< OneD, const NekDouble > &  soln 
)

Calculate the $ L^2 $ error of the $ Q_{\rm dir} $ derivative using the consistent DG evaluation of $ Q_{\rm dir} $.

The solution provided is of the primative variation at the quadrature points and the derivative is compared to the discrete derivative at these points, which is likely to be undesirable unless using a much higher number of quadrature points than the polynomial order used to evaluate $ Q_{\rm dir} $.

Definition at line 1680 of file DisContField2D.cpp.

References Nektar::MultiRegions::ExpList::BwdTrans(), Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::MultiRegions::ExpList::L2(), Nektar::MultiRegions::ExpList::m_coeff_offset, Nektar::MultiRegions::ExpList::m_coeffs, Nektar::MultiRegions::ExpList::m_ncoeffs, Nektar::MultiRegions::ExpList::m_npoints, Nektar::MultiRegions::ExpList::m_offset_elmt_id, Nektar::MultiRegions::ExpList::m_phys, m_trace, m_traceMap, Vmath::Vcopy(), and Vmath::Vsub().

{
int i,e,ncoeff_edge;
Array<OneD, const NekDouble> tmp_coeffs;
Array<OneD, NekDouble> out_d(m_ncoeffs), out_tmp;
Array<OneD, Array< OneD, StdRegions::StdExpansionSharedPtr> >
&elmtToTrace = m_traceMap->GetElmtToTrace();
int eid,cnt;
int LocBndCoeffs = m_traceMap->GetNumLocalBndCoeffs();
Array<OneD, NekDouble> loc_lambda(LocBndCoeffs), edge_lambda;
m_traceMap->GlobalToLocalBnd(m_trace->GetCoeffs(),loc_lambda);
edge_lambda = loc_lambda;
// Calculate Q using standard DG formulation.
for(i = cnt = 0; i < GetExpSize(); ++i)
{
eid = m_offset_elmt_id[i];
// Probably a better way of setting up lambda than this.
// Note cannot use PutCoeffsInToElmts since lambda space
// is mapped during the solve.
int nEdges = (*m_exp)[i]->GetNedges();
Array<OneD, Array<OneD, NekDouble> > edgeCoeffs(nEdges);
for(e = 0; e < nEdges; ++e)
{
edgedir = (*m_exp)[eid]->GetEorient(e);
ncoeff_edge = elmtToTrace[eid][e]->GetNcoeffs();
edgeCoeffs[e] = Array<OneD, NekDouble>(ncoeff_edge);
Vmath::Vcopy(ncoeff_edge, edge_lambda, 1, edgeCoeffs[e], 1);
elmtToTrace[eid][e]->SetCoeffsToOrientation(
edgedir, edgeCoeffs[e], edgeCoeffs[e]);
edge_lambda = edge_lambda + ncoeff_edge;
}
(*m_exp)[eid]->DGDeriv(dir,
tmp_coeffs=m_coeffs+m_coeff_offset[eid],
elmtToTrace[eid],
edgeCoeffs,
out_tmp = out_d+cnt);
cnt += (*m_exp)[eid]->GetNcoeffs();
}
BwdTrans(out_d,m_phys);
return L2(m_phys);
}
bool Nektar::MultiRegions::DisContField2D::SameTypeOfBoundaryConditions ( const DisContField2D In)
protected

For each boundary region, checks that the types and number of boundary expansions in that region match.

Parameters
InContField2D to compare with.
Returns
True if boundary conditions match.

Definition at line 559 of file DisContField2D.cpp.

References Nektar::MultiRegions::ExpList::GetExpSize(), m_bndCondExpansions, m_bndConditions, Nektar::MultiRegions::ExpList::m_comm, and Nektar::LibUtilities::ReduceMin.

Referenced by Nektar::MultiRegions::ContField2D::ContField2D(), and DisContField2D().

{
int i;
bool returnval = true;
for(i = 0; i < m_bndConditions.num_elements(); ++i)
{
// check to see if boundary condition type is the same
// and there are the same number of boundary
// conditions in the boundary definition.
if((m_bndConditions[i]->GetBoundaryConditionType()
!= In.m_bndConditions[i]->GetBoundaryConditionType())||
!= In.m_bndCondExpansions[i]->GetExpSize()))
{
returnval = false;
break;
}
}
// Compare with all other processes. Return true only if all
// processes report having the same boundary conditions.
int vSame = (returnval?1:0);
m_comm->GetRowComm()->AllReduce(vSame, LibUtilities::ReduceMin);
return (vSame == 1);
}
void Nektar::MultiRegions::DisContField2D::SetUpDG ( const std::string  variable = "DefaultVar")
protected

Set up all DG member variables and maps.

Definition at line 355 of file DisContField2D.cpp.

References Nektar::StdRegions::StdExpansion::as(), ASSERTL1, Nektar::StdRegions::eForwards, Nektar::SpatialDomains::ePeriodic, Nektar::MultiRegions::_PeriodicEntity::id, IsLeftAdjacentEdge(), Nektar::MultiRegions::_PeriodicEntity::isLocal, Nektar::iterator, m_bndCondExpansions, m_bndConditions, m_boundaryEdges, Nektar::MultiRegions::ExpList::m_exp, m_globalBndMat, Nektar::MultiRegions::ExpList::m_graph, m_leftAdjacentEdges, m_periodicBwdCopy, m_periodicEdges, m_periodicFwdCopy, Nektar::MultiRegions::ExpList::m_session, m_trace, m_traceMap, Nektar::MultiRegions::NullExpListSharedPtr, Nektar::MultiRegions::_PeriodicEntity::orient, and Nektar::MultiRegions::ExpList::SetUpPhysNormals().

Referenced by DisContField2D(), and v_GetTrace().

{
// Check for multiple calls
{
return;
}
boost::dynamic_pointer_cast<SpatialDomains::MeshGraph2D>(
// Set up matrix map
// Set up trace space
graph2D, m_periodicEdges);
m_trace = boost::dynamic_pointer_cast<ExpList>(trace);
AllocateSharedPtr(m_session, graph2D, trace, *this,
variable);
Array<OneD, Array<OneD, StdRegions::StdExpansionSharedPtr> >
&elmtToTrace = m_traceMap->GetElmtToTrace();
// Scatter trace segments to 2D elements. For each element, we find
// the trace segment associated to each edge. The element then
// retains a pointer to the trace space segments, to ensure
// uniqueness of normals when retrieving from two adjoining elements
// which do not lie in a plane.
for (int i = 0; i < m_exp->size(); ++i)
{
for (int j = 0; j < (*m_exp)[i]->GetNedges(); ++j)
{
(*m_exp)[i]->as<LocalRegions::Expansion2D>();
elmtToTrace[i][j]->as<LocalRegions::Expansion1D>();
elmtToTrace[i][j]->as<LocalRegions::Expansion> ();
exp2d->SetEdgeExp (j, exp );
exp1d->SetAdjacentElementExp(j, exp2d);
}
}
// Set up physical normals
// Set up information for parallel and periodic problems.
for (int i = 0; i < m_trace->GetExpSize(); ++i)
{
int offset = m_trace->GetPhys_Offset(i);
int traceGeomId = traceEl->GetGeom1D()->GetGlobalID();
traceGeomId);
if (pIt != m_periodicEdges.end() && !pIt->second[0].isLocal)
{
if (traceGeomId != min(pIt->second[0].id, traceGeomId))
{
traceEl->GetLeftAdjacentElementExp()->NegateEdgeNormal(
traceEl->GetLeftAdjacentElementEdge());
}
}
else if (m_traceMap->GetTraceToUniversalMapUnique(offset) < 0)
{
traceEl->GetLeftAdjacentElementExp()->NegateEdgeNormal(
traceEl->GetLeftAdjacentElementEdge());
}
}
int cnt, n, e;
// Identify boundary edges
for(cnt = 0, n = 0; n < m_bndCondExpansions.num_elements(); ++n)
{
if (m_bndConditions[n]->GetBoundaryConditionType() !=
{
for(e = 0; e < m_bndCondExpansions[n]->GetExpSize(); ++e)
{
m_traceMap->GetBndCondTraceToGlobalTraceMap(cnt+e));
}
}
cnt += m_bndCondExpansions[n]->GetExpSize();
}
// Set up information for periodic boundary conditions.
boost::unordered_map<int,pair<int,int> > perEdgeToExpMap;
boost::unordered_map<int,pair<int,int> >::iterator it2;
for (cnt = n = 0; n < m_exp->size(); ++n)
{
for (e = 0; e < (*m_exp)[n]->GetNedges(); ++e, ++cnt)
{
(*m_exp)[n]->GetGeom()->GetEid(e));
if (it != m_periodicEdges.end())
{
perEdgeToExpMap[it->first] = make_pair(n, e);
}
}
}
// Set up left-adjacent edge list.
m_leftAdjacentEdges.resize(cnt);
cnt = 0;
for (int i = 0; i < m_exp->size(); ++i)
{
for (int j = 0; j < (*m_exp)[i]->GetNedges(); ++j, ++cnt)
{
}
}
// Set up mapping to copy Fwd of periodic bcs to Bwd of other edge.
cnt = 0;
for (int n = 0; n < m_exp->size(); ++n)
{
for (int e = 0; e < (*m_exp)[n]->GetNedges(); ++e, ++cnt)
{
int edgeGeomId = (*m_exp)[n]->GetGeom()->GetEid(e);
int offset = m_trace->GetPhys_Offset(
elmtToTrace[n][e]->GetElmtId());
// Check to see if this face is periodic.
PeriodicMap::iterator it = m_periodicEdges.find(edgeGeomId);
if (it != m_periodicEdges.end())
{
const PeriodicEntity &ent = it->second[0];
it2 = perEdgeToExpMap.find(ent.id);
if (it2 == perEdgeToExpMap.end())
{
if (m_session->GetComm()->GetRowComm()->GetSize() > 1 &&
!ent.isLocal)
{
continue;
}
else
{
ASSERTL1(false, "Periodic edge not found!");
}
}
"Periodic edge in non-forward space?");
int offset2 = m_trace->GetPhys_Offset(
elmtToTrace[it2->second.first][it2->second.second]->
GetElmtId());
// Calculate relative orientations between edges to
// calculate copying map.
int nquad = elmtToTrace[n][e]->GetNumPoints(0);
vector<int> tmpBwd(nquad);
vector<int> tmpFwd(nquad);
if (ent.orient == StdRegions::eForwards)
{
for (int i = 0; i < nquad; ++i)
{
tmpBwd[i] = offset2 + i;
tmpFwd[i] = offset + i;
}
}
else
{
for (int i = 0; i < nquad; ++i)
{
tmpBwd[i] = offset2 + i;
tmpFwd[i] = offset + nquad - i - 1;
}
}
for (int i = 0; i < nquad; ++i)
{
m_periodicFwdCopy.push_back(tmpFwd[i]);
m_periodicBwdCopy.push_back(tmpBwd[i]);
}
}
}
}
}
void Nektar::MultiRegions::DisContField2D::v_AddFwdBwdTraceIntegral ( const Array< OneD, const NekDouble > &  Fwd,
const Array< OneD, const NekDouble > &  Bwd,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Add trace contributions into elemental coefficient spaces.

Given some quantity $ \vec{q} $, calculate the elemental integral

\[ \int_{\Omega^e} \vec{q}, \mathrm{d}S \]

and adds this to the coefficient space provided by outarray. The value of q is determined from the routine IsLeftAdjacentEdge() which if true we use Fwd else we use Bwd

See Also
Expansion2D::AddEdgeNormBoundaryInt
Parameters
FwdThe trace quantities associated with left (fwd) adjancent elmt.
BwdThe trace quantities associated with right (bwd) adjacent elet.
outarrayResulting 2D coefficient space.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1574 of file DisContField2D.cpp.

References Nektar::MultiRegions::ExpList::GetCoeff_Offset(), Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::MultiRegions::ExpList::GetTrace(), IsLeftAdjacentEdge(), and m_traceMap.

{
int e,n,offset, t_offset;
Array<OneD, NekDouble> e_outarray;
Array<OneD, Array<OneD, StdRegions::StdExpansionSharedPtr> >
&elmtToTrace = m_traceMap->GetElmtToTrace();
for (n = 0; n < GetExpSize(); ++n)
{
offset = GetCoeff_Offset(n);
for (e = 0; e < (*m_exp)[n]->GetNedges(); ++e)
{
t_offset = GetTrace()->GetPhys_Offset(
elmtToTrace[n][e]->GetElmtId());
// Evaluate upwind flux less local edge
if (IsLeftAdjacentEdge(n, e))
{
(*m_exp)[n]->AddEdgeNormBoundaryInt(
e, elmtToTrace[n][e], Fwd+t_offset,
e_outarray = outarray+offset);
}
else
{
(*m_exp)[n]->AddEdgeNormBoundaryInt(
e, elmtToTrace[n][e], Bwd+t_offset,
e_outarray = outarray+offset);
}
}
}
}
void Nektar::MultiRegions::DisContField2D::v_AddTraceIntegral ( const Array< OneD, const NekDouble > &  Fx,
const Array< OneD, const NekDouble > &  Fy,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1484 of file DisContField2D.cpp.

References Nektar::MultiRegions::ExpList::GetCoeff_Offset(), Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::MultiRegions::ExpList::GetTrace(), and m_traceMap.

{
int e, n, offset, t_offset;
Array<OneD, NekDouble> e_outarray;
Array<OneD, Array<OneD, StdRegions::StdExpansionSharedPtr> >
&elmtToTrace = m_traceMap->GetElmtToTrace();
for(n = 0; n < GetExpSize(); ++n)
{
offset = GetCoeff_Offset(n);
for(e = 0; e < (*m_exp)[n]->GetNedges(); ++e)
{
t_offset = GetTrace()->GetPhys_Offset(
elmtToTrace[n][e]->GetElmtId());
(*m_exp)[n]->AddEdgeNormBoundaryInt(
e,elmtToTrace[n][e],
Fx + t_offset,
Fy + t_offset,
e_outarray = outarray+offset);
}
}
}
void Nektar::MultiRegions::DisContField2D::v_AddTraceIntegral ( const Array< OneD, const NekDouble > &  Fn,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Add trace contributions into elemental coefficient spaces.

Given some quantity $ \vec{Fn} $, which conatins this routine calculates the integral

\[ \int_{\Omega^e} \vec{Fn}, \mathrm{d}S \]

and adds this to the coefficient space provided by outarray.

See Also
Expansion2D::AddEdgeNormBoundaryInt
Parameters
FnThe trace quantities.
outarrayResulting 2D coefficient space.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1528 of file DisContField2D.cpp.

References Nektar::MultiRegions::ExpList::GetCoeff_Offset(), Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::MultiRegions::ExpList::GetTrace(), and m_traceMap.

{
int e, n, offset, t_offset;
Array<OneD, NekDouble> e_outarray;
Array<OneD, Array<OneD, StdRegions::StdExpansionSharedPtr> >
&elmtToTrace = m_traceMap->GetElmtToTrace();
for(n = 0; n < GetExpSize(); ++n)
{
offset = GetCoeff_Offset(n);
for(e = 0; e < (*m_exp)[n]->GetNedges(); ++e)
{
t_offset = GetTrace()->GetPhys_Offset(
elmtToTrace[n][e]->GetElmtId());
(*m_exp)[n]->AddEdgeNormBoundaryInt(
e, elmtToTrace[n][e], Fn+t_offset,
e_outarray = outarray+offset);
}
}
}
void Nektar::MultiRegions::DisContField2D::v_EvaluateBoundaryConditions ( const NekDouble  time = 0.0,
const std::string  varName = "",
const NekDouble  x2_in = NekConstants::kNekUnsetDouble,
const NekDouble  x3_in = NekConstants::kNekUnsetDouble 
)
protectedvirtual

Evaluates the boundary condition expansions, bndCondExpansions, given the information provided by bndConditions.

Parameters
timeThe time at which the boundary conditions should be evaluated.
bndCondExpansionsList of boundary conditions.
bndConditionsInformation about the boundary conditions.

This will only be undertaken for time dependent boundary conditions unless time == 0.0 which is the case when the method is called from the constructor.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 2115 of file DisContField2D.cpp.

References ASSERTL0, Nektar::SpatialDomains::eDirichlet, Nektar::SpatialDomains::eNeumann, Nektar::SpatialDomains::eRobin, Nektar::SpatialDomains::eTimeDependent, Nektar::LibUtilities::Equation::Evaluate(), Nektar::MultiRegions::ExpList::ExtractFileBCs(), Vmath::Fill(), Nektar::NekConstants::kNekUnsetDouble, m_bndCondExpansions, and m_bndConditions.

{
int i;
int npoints;
int nbnd = m_bndCondExpansions.num_elements();
for (i = 0; i < nbnd; ++i)
{
if (time == 0.0 ||
m_bndConditions[i]->GetUserDefined() ==
{
locExpList = m_bndCondExpansions[i];
npoints = locExpList->GetNpoints();
Array<OneD, NekDouble> x0(npoints, 0.0);
Array<OneD, NekDouble> x1(npoints, 0.0);
Array<OneD, NekDouble> x2(npoints, 0.0);
// Homogeneous input case for x2.
{
locExpList->GetCoords(x0, x1, x2);
}
else
{
locExpList->GetCoords(x0, x1, x2);
Vmath::Fill(npoints, x2_in, x2, 1);
}
if (m_bndConditions[i]->GetBoundaryConditionType()
{
string filebcs = boost::static_pointer_cast<
m_bndConditions[i])->m_filename;
if (filebcs != "")
{
ExtractFileBCs(filebcs, varName, locExpList);
}
else
{
boost::static_pointer_cast<
m_dirichletCondition;
condition.Evaluate(x0, x1, x2, time,
locExpList->UpdatePhys());
}
locExpList->FwdTrans_BndConstrained(
locExpList->GetPhys(),
locExpList->UpdateCoeffs());
}
else if (m_bndConditions[i]->GetBoundaryConditionType()
{
string filebcs = boost::static_pointer_cast<
m_bndConditions[i])->m_filename;
if (filebcs != "")
{
ExtractFileBCs(filebcs, varName, locExpList);
}
else
{
boost::static_pointer_cast<
m_neumannCondition;
condition.Evaluate(x0, x1, x2, time,
locExpList->UpdatePhys());
}
locExpList->IProductWRTBase(
locExpList->GetPhys(),
locExpList->UpdateCoeffs());
}
else if (m_bndConditions[i]->GetBoundaryConditionType()
{
string filebcs = boost::static_pointer_cast<
(m_bndConditions[i])->m_filename;
if (filebcs != "")
{
ExtractFileBCs(filebcs, varName, locExpList);
}
else
{
boost::static_pointer_cast<
m_robinFunction;
condition.Evaluate(x0, x1, x2, time,
locExpList->UpdatePhys());
}
boost::static_pointer_cast<
m_bndConditions[i])->m_robinPrimitiveCoeff;
locExpList->IProductWRTBase(
locExpList->GetPhys(),
locExpList->UpdateCoeffs());
// put primitive coefficient into the physical
// space storage
coeff.Evaluate(x0, x1, x2, time,
locExpList->UpdatePhys());
}
else
{
ASSERTL0(false, "This type of BC not implemented yet");
}
}
}
}
void Nektar::MultiRegions::DisContField2D::v_ExtractTracePhys ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

This method extracts the trace (edges in 2D) from the field inarray and puts the values in outarray.

It assumes the field is C0 continuous so that it can overwrite the edge data when visited by the two adjacent elements.

Parameters
inarrayAn array containing the 2D data from which we wish to extract the edge data.
outarrayThe resulting edge information.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1454 of file DisContField2D.cpp.

References ASSERTL1, Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::MultiRegions::ExpList::GetPhys_Offset(), m_trace, and m_traceMap.

Referenced by v_ExtractTracePhys().

{
// Loop over elemente and collect forward expansion
int nexp = GetExpSize();
int n, e, offset, phys_offset;
Array<OneD,NekDouble> e_tmp;
Array<OneD, Array<OneD, StdRegions::StdExpansionSharedPtr> >
&elmtToTrace = m_traceMap->GetElmtToTrace();
ASSERTL1(outarray.num_elements() >= m_trace->GetNpoints(),
"input array is of insufficient length");
// use m_trace tmp space in element to fill values
for(n = 0; n < nexp; ++n)
{
phys_offset = GetPhys_Offset(n);
for(e = 0; e < (*m_exp)[n]->GetNedges(); ++e)
{
offset = m_trace->GetPhys_Offset(
elmtToTrace[n][e]->GetElmtId());
(*m_exp)[n]->GetEdgePhysVals(e, elmtToTrace[n][e],
inarray + phys_offset,
e_tmp = outarray + offset);
}
}
}
void Nektar::MultiRegions::DisContField2D::v_ExtractTracePhys ( Array< OneD, NekDouble > &  outarray)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1434 of file DisContField2D.cpp.

References ASSERTL1, Nektar::MultiRegions::ExpList::m_phys, Nektar::MultiRegions::ExpList::m_physState, and v_ExtractTracePhys().

{
"Field must be in physical state to extract trace space.");
}
void Nektar::MultiRegions::DisContField2D::v_FillBndCondFromField ( void  )
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1429 of file DisContField2D.cpp.

{
}
void Nektar::MultiRegions::DisContField2D::v_GeneralMatrixOp ( const GlobalMatrixKey gkey,
const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
CoeffState  coeffstate = eLocal 
)
protectedvirtual

Calculates the result of the multiplication of a global matrix of type specified by mkey with a vector given by inarray.

Parameters
mkeyKey representing desired matrix multiplication.
inarrayInput vector.
outarrayResulting multiplication.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1876 of file DisContField2D.cpp.

References Nektar::eWrapper, Nektar::MultiRegions::ExpList::GetBlockMatrix(), and m_traceMap.

{
int LocBndCoeffs = m_traceMap->GetNumLocalBndCoeffs();
Array<OneD, NekDouble> loc_lambda(LocBndCoeffs);
DNekVec LocLambda(LocBndCoeffs,loc_lambda,eWrapper);
const DNekScalBlkMatSharedPtr& HDGHelm = GetBlockMatrix(gkey);
m_traceMap->GlobalToLocalBnd(inarray, loc_lambda);
LocLambda = (*HDGHelm) * LocLambda;
m_traceMap->AssembleBnd(loc_lambda,outarray);
}
virtual const Array<OneD,const MultiRegions::ExpListSharedPtr>& Nektar::MultiRegions::DisContField2D::v_GetBndCondExpansions ( void  )
inlineprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 227 of file DisContField2D.h.

References m_bndCondExpansions.

{
}
virtual const Array<OneD,const SpatialDomains::BoundaryConditionShPtr>& Nektar::MultiRegions::DisContField2D::v_GetBndConditions ( void  )
inlineprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 234 of file DisContField2D.h.

References m_bndConditions.

{
}
void Nektar::MultiRegions::DisContField2D::v_GetBoundaryToElmtMap ( Array< OneD, int > &  ElmtID,
Array< OneD, int > &  EdgeID 
)
protectedvirtual

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

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1614 of file DisContField2D.cpp.

References Nektar::MultiRegions::ExpList::GetExpSize(), m_bndCondExpansions, m_bndConditions, and Nektar::MultiRegions::ExpList::m_graph.

{
map<int, int> globalIdMap;
int i,n;
int cnt;
int nbcs = 0;
boost::dynamic_pointer_cast<SpatialDomains::MeshGraph2D>(
// Populate global ID map (takes global geometry ID to local
// expansion list ID).
for (i = 0; i < GetExpSize(); ++i)
{
globalIdMap[(*m_exp)[i]->GetGeom()->GetGlobalID()] = i;
}
// Determine number of boundary condition expansions.
for(i = 0; i < m_bndConditions.num_elements(); ++i)
{
nbcs += m_bndCondExpansions[i]->GetExpSize();
}
// make sure arrays are of sufficient length
if(ElmtID.num_elements() != nbcs)
{
ElmtID = Array<OneD, int>(nbcs);
}
if(EdgeID.num_elements() != nbcs)
{
EdgeID = Array<OneD, int>(nbcs);
}
for (cnt = n = 0; n < m_bndCondExpansions.num_elements(); ++n)
{
for (i = 0; i < m_bndCondExpansions[n]->GetExpSize();
++i, ++cnt)
{
exp1d = m_bndCondExpansions[n]->GetExp(i)->
as<LocalRegions::Expansion1D>();
// Use edge to element map from MeshGraph2D.
graph2D->GetElementsFromEdge(exp1d->GetGeom1D());
ElmtID[cnt] = globalIdMap[(*tmp)[0]->
m_Element->GetGlobalID()];
EdgeID[cnt] = (*tmp)[0]->m_EdgeIndx;
}
}
}
void Nektar::MultiRegions::DisContField2D::v_GetFwdBwdTracePhys ( const Array< OneD, const NekDouble > &  field,
Array< OneD, NekDouble > &  Fwd,
Array< OneD, NekDouble > &  Bwd 
)
protectedvirtual

This method extracts the "forward" and "backward" trace data from the array field and puts the data into output vectors Fwd and Bwd.

We first define the convention which defines "forwards" and "backwards". First an association is made between the edge of each element and its corresponding edge in the trace space using the mapping m_traceMap. The element can either be left-adjacent or right-adjacent to this trace edge (see Expansion1D::GetLeftAdjacentElementExp). Boundary edges are always left-adjacent since left-adjacency is populated first.

If the element is left-adjacent we extract the edge trace data from field into the forward trace space Fwd; otherwise, we place it in the backwards trace space Bwd. In this way, we form a unique set of trace normals since these are always extracted from left-adjacent elements.

Parameters
fieldis a NekDouble array which contains the 2D data from which we wish to extract the backward and forward orientated trace/edge arrays.
FwdThe resulting forwards space.
BwdThe resulting backwards space.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1325 of file DisContField2D.cpp.

References ASSERTL0, Nektar::SpatialDomains::eDirichlet, Nektar::SpatialDomains::eNeumann, Nektar::SpatialDomains::ePeriodic, Nektar::SpatialDomains::eRobin, Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::MultiRegions::ExpList::GetPhys(), Nektar::MultiRegions::ExpList::GetPhys_Offset(), Nektar::iterator, m_bndCondExpansions, m_bndConditions, m_leftAdjacentEdges, m_periodicBwdCopy, m_periodicFwdCopy, m_trace, m_traceMap, npts, Vmath::Vcopy(), and Vmath::Zero().

Referenced by v_GetFwdBwdTracePhys().

{
// Loop over elements and collect forward expansion
int nexp = GetExpSize();
int cnt, n, e, npts, phys_offset;
Array<OneD,NekDouble> e_tmp;
boost::unordered_map<int,pair<int,int> >::iterator it3;
Array<OneD, Array<OneD, StdRegions::StdExpansionSharedPtr> >
&elmtToTrace = m_traceMap->GetElmtToTrace();
// Zero forward/backward vectors.
Vmath::Zero(Fwd.num_elements(), Fwd, 1);
Vmath::Zero(Bwd.num_elements(), Bwd, 1);
for(cnt = n = 0; n < nexp; ++n)
{
exp2d = (*m_exp)[n]->as<LocalRegions::Expansion2D>();
phys_offset = GetPhys_Offset(n);
for(e = 0; e < exp2d->GetNedges(); ++e, ++cnt)
{
int offset = m_trace->GetPhys_Offset(
elmtToTrace[n][e]->GetElmtId());
{
exp2d->GetEdgePhysVals(e, elmtToTrace[n][e],
field + phys_offset,
e_tmp = Fwd + offset);
}
else
{
exp2d->GetEdgePhysVals(e, elmtToTrace[n][e],
field + phys_offset,
e_tmp = Bwd + offset);
}
}
}
// Fill boundary conditions into missing elements.
int id1, id2 = 0;
for(cnt = n = 0; n < m_bndCondExpansions.num_elements(); ++n)
{
if (m_bndConditions[n]->GetBoundaryConditionType() ==
{
for(e = 0; e < m_bndCondExpansions[n]->GetExpSize(); ++e)
{
npts = m_bndCondExpansions[n]->GetExp(e)->GetNumPoints(0);
id1 = m_bndCondExpansions[n]->GetPhys_Offset(e);
id2 = m_trace->GetPhys_Offset(
m_traceMap->GetBndCondTraceToGlobalTraceMap(cnt+e));
&(m_bndCondExpansions[n]->GetPhys())[id1], 1,
&Bwd[id2], 1);
}
cnt += e;
}
else if (m_bndConditions[n]->GetBoundaryConditionType() ==
m_bndConditions[n]->GetBoundaryConditionType() ==
{
for(e = 0; e < m_bndCondExpansions[n]->GetExpSize(); ++e)
{
npts = m_bndCondExpansions[n]->GetExp(e)->GetNumPoints(0);
id1 = m_bndCondExpansions[n]->GetPhys_Offset(e);
"Method not set up for non-zero Neumann "
"boundary condition");
id2 = m_trace->GetPhys_Offset(
m_traceMap->GetBndCondTraceToGlobalTraceMap(cnt+e));
Vmath::Vcopy(npts, &Fwd[id2], 1, &Bwd[id2], 1);
}
cnt += e;
}
else if (m_bndConditions[n]->GetBoundaryConditionType() !=
{
ASSERTL0(false,
"Method not set up for this boundary condition.");
}
}
// Copy any periodic boundary conditions.
for (n = 0; n < m_periodicFwdCopy.size(); ++n)
{
}
// Do parallel exchange for forwards/backwards spaces.
m_traceMap->UniversalTraceAssemble(Fwd);
m_traceMap->UniversalTraceAssemble(Bwd);
}
void Nektar::MultiRegions::DisContField2D::v_GetFwdBwdTracePhys ( Array< OneD, NekDouble > &  Fwd,
Array< OneD, NekDouble > &  Bwd 
)
protectedvirtual
virtual void Nektar::MultiRegions::DisContField2D::v_GetPeriodicEntities ( PeriodicMap periodicVerts,
PeriodicMap periodicEdges,
PeriodicMap periodicFaces 
)
inlineprotectedvirtual

Obtain a copy of the periodic edges and vertices for this field.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 211 of file DisContField2D.h.

References m_periodicEdges, and m_periodicVerts.

{
periodicVerts = m_periodicVerts;
periodicEdges = m_periodicEdges;
}
map< int, RobinBCInfoSharedPtr > Nektar::MultiRegions::DisContField2D::v_GetRobinBCInfo ( void  )
protectedvirtual

Search through the edge expansions and identify which ones have Robin/Mixed type boundary conditions.

If a Robin boundary is found then store the edge ID of the boundary condition and the array of points of the physical space boundary condition which are hold the boundary condition primitive variable coefficient at the quatrature points

Returns
A map containing the Robin boundary condition information using a key of the element ID.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1904 of file DisContField2D.cpp.

References Nektar::SpatialDomains::eRobin, Nektar::MultiRegions::ExpList::GetBoundaryToElmtMap(), m_bndCondExpansions, and m_bndConditions.

{
int i,cnt;
map<int, RobinBCInfoSharedPtr> returnval;
Array<OneD, int> ElmtID,EdgeID;
GetBoundaryToElmtMap(ElmtID,EdgeID);
for(cnt = i = 0; i < m_bndCondExpansions.num_elements(); ++i)
{
if(m_bndConditions[i]->GetBoundaryConditionType() ==
{
int e,elmtid;
Array<OneD, NekDouble> Array_tmp;
locExpList = m_bndCondExpansions[i];
for(e = 0; e < locExpList->GetExpSize(); ++e)
{
EdgeID[cnt+e],
Array_tmp = locExpList->GetPhys() +
locExpList->GetPhys_Offset(e));
elmtid = ElmtID[cnt+e];
// make link list if necessary
if(returnval.count(elmtid) != 0)
{
rInfo->next = returnval.find(elmtid)->second;
}
returnval[elmtid] = rInfo;
}
}
cnt += m_bndCondExpansions[i]->GetExpSize();
}
return returnval;
}
virtual ExpListSharedPtr& Nektar::MultiRegions::DisContField2D::v_GetTrace ( )
inlinevirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 86 of file DisContField2D.h.

References m_trace, Nektar::MultiRegions::NullExpListSharedPtr, and SetUpDG().

{
{
}
return m_trace;
}
virtual AssemblyMapDGSharedPtr& Nektar::MultiRegions::DisContField2D::v_GetTraceMap ( void  )
inlineprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 221 of file DisContField2D.h.

References m_traceMap.

{
return m_traceMap;
}
void Nektar::MultiRegions::DisContField2D::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 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1735 of file DisContField2D.cpp.

References Nektar::SpatialDomains::eDirichlet, Nektar::StdRegions::eHybridDGHelmBndLam, Nektar::StdRegions::eHybridDGLamToU, Nektar::StdRegions::eInvHybridDGHelmholtz, Nektar::eWrapper, Nektar::MultiRegions::ExpList::GetBlockMatrix(), Nektar::MultiRegions::ExpList::GetExpSize(), GetGlobalBndLinSys(), Nektar::MultiRegions::ExpList::GetNcoeffs(), Nektar::MultiRegions::ExpList::IProductWRTBase(), m_bndCondExpansions, m_bndConditions, Nektar::MultiRegions::ExpList::m_ncoeffs, Nektar::MultiRegions::ExpList::m_offset_elmt_id, m_trace, m_traceMap, Vmath::Neg(), Nektar::MultiRegions::NullAssemblyMapSharedPtr, Nektar::Transpose(), and Vmath::Zero().

{
int i,j,n,cnt,cnt1,nbndry;
int nexp = GetExpSize();
Array<OneD,NekDouble> f(m_ncoeffs);
Array<OneD,NekDouble> e_f, e_l;
//----------------------------------
// Setup RHS Inner product
//----------------------------------
IProductWRTBase(inarray,f);
//----------------------------------
// Solve continuous flux System
//----------------------------------
int GloBndDofs = m_traceMap->GetNumGlobalBndCoeffs();
int NumDirichlet = m_traceMap->GetNumLocalDirBndCoeffs();
int e_ncoeffs,id;
// Retrieve block matrix of U^e
GlobalMatrixKey HDGLamToUKey(StdRegions::eHybridDGLamToU,
NullAssemblyMapSharedPtr,factors,varcoeff);
HDGLamToUKey);
// Retrieve global trace space storage, \Lambda, from trace
// expansion
Array<OneD,NekDouble> BndSol = m_trace->UpdateCoeffs();
// Create trace space forcing, F
Array<OneD,NekDouble> BndRhs(GloBndDofs,0.0);
// Zero \Lambda
Vmath::Zero(GloBndDofs,BndSol,1);
// Retrieve number of local trace space coefficients N_{\lambda},
// and set up local elemental trace solution \lambda^e.
int LocBndCoeffs = m_traceMap->GetNumLocalBndCoeffs();
Array<OneD, NekDouble> loc_lambda(LocBndCoeffs);
DNekVec LocLambda(LocBndCoeffs,loc_lambda,eWrapper);
//----------------------------------
// Evaluate Trace Forcing vector F
// Kirby et al, 2010, P23, Step 5.
//----------------------------------
// Loop over all expansions in the domain
for(cnt = cnt1 = n = 0; n < nexp; ++n)
{
nbndry = (*m_exp)[m_offset_elmt_id[n]]->NumDGBndryCoeffs();
e_ncoeffs = (*m_exp)[m_offset_elmt_id[n]]->GetNcoeffs();
e_f = f + cnt;
e_l = loc_lambda + cnt1;
// Local trace space \lambda^e
DNekVec Floc (nbndry, e_l, eWrapper);
// Local forcing f^e
DNekVec ElmtFce (e_ncoeffs, e_f, eWrapper);
// Compute local (U^e)^{\top} f^e
Floc = Transpose(*(HDGLamToU->GetBlock(n,n)))*ElmtFce;
cnt += e_ncoeffs;
cnt1 += nbndry;
}
// Assemble local \lambda_e into global \Lambda
m_traceMap->AssembleBnd(loc_lambda,BndRhs);
// Copy Dirichlet boundary conditions and weak forcing into trace
// space
cnt = 0;
for(i = 0; i < m_bndCondExpansions.num_elements(); ++i)
{
if(m_bndConditions[i]->GetBoundaryConditionType() ==
{
for(j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); ++j)
{
id = m_traceMap->GetBndCondCoeffsToGlobalCoeffsMap(cnt++);
BndSol[id] = m_bndCondExpansions[i]->GetCoeffs()[j];
}
}
else
{
//Add weak boundary condition to trace forcing
for(j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); ++j)
{
id = m_traceMap->GetBndCondCoeffsToGlobalCoeffsMap(cnt++);
BndRhs[id] += m_bndCondExpansions[i]->GetCoeffs()[j];
}
}
}
//----------------------------------
// Solve trace problem: \Lambda = K^{-1} F
// K is the HybridDGHelmBndLam matrix.
//----------------------------------
if(GloBndDofs - NumDirichlet > 0)
{
GlobalLinSysKey key(StdRegions::eHybridDGHelmBndLam,
m_traceMap,factors,varcoeff);
LinSys->Solve(BndRhs,BndSol,m_traceMap);
}
//----------------------------------
// Internal element solves
//----------------------------------
GlobalMatrixKey invHDGhelmkey(StdRegions::eInvHybridDGHelmholtz,
NullAssemblyMapSharedPtr,factors,varcoeff);
invHDGhelmkey);
DNekVec out(m_ncoeffs,outarray,eWrapper);
Vmath::Zero(m_ncoeffs,outarray,1);
// get local trace solution from BndSol
m_traceMap->GlobalToLocalBnd(BndSol,loc_lambda);
// out = u_f + u_lam = (*InvHDGHelm)*f + (LamtoU)*Lam
out = (*InvHDGHelm)*F + (*HDGLamToU)*LocLambda;
}
virtual MultiRegions::ExpListSharedPtr& Nektar::MultiRegions::DisContField2D::v_UpdateBndCondExpansion ( int  i)
inlineprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 240 of file DisContField2D.h.

References m_bndCondExpansions.

{
}
virtual Array<OneD, SpatialDomains::BoundaryConditionShPtr>& Nektar::MultiRegions::DisContField2D::v_UpdateBndConditions ( )
inlineprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 246 of file DisContField2D.h.

References m_bndConditions.

{
}

Member Data Documentation

Array<OneD,MultiRegions::ExpListSharedPtr> Nektar::MultiRegions::DisContField2D::m_bndCondExpansions
protected
Array<OneD,SpatialDomains::BoundaryConditionShPtr> Nektar::MultiRegions::DisContField2D::m_bndConditions
protected
std::set<int> Nektar::MultiRegions::DisContField2D::m_boundaryEdges
protected

A set storing the global IDs of any boundary edges.

Definition at line 127 of file DisContField2D.h.

Referenced by DisContField2D(), IsLeftAdjacentEdge(), and SetUpDG().

Array<OneD,StdRegions::Orientation> Nektar::MultiRegions::DisContField2D::m_edgedir
protected

Definition at line 122 of file DisContField2D.h.

GlobalLinSysMapShPtr Nektar::MultiRegions::DisContField2D::m_globalBndMat
protected

Definition at line 116 of file DisContField2D.h.

Referenced by DisContField2D(), GetGlobalBndLinSys(), and SetUpDG().

vector<bool> Nektar::MultiRegions::DisContField2D::m_leftAdjacentEdges
protected

Definition at line 152 of file DisContField2D.h.

Referenced by DisContField2D(), SetUpDG(), and v_GetFwdBwdTracePhys().

Array<OneD, Array<OneD, unsigned int> > Nektar::MultiRegions::DisContField2D::m_mapEdgeToElmn
protected

Definition at line 120 of file DisContField2D.h.

vector<int> Nektar::MultiRegions::DisContField2D::m_periodicBwdCopy
protected

Definition at line 146 of file DisContField2D.h.

Referenced by DisContField2D(), SetUpDG(), and v_GetFwdBwdTracePhys().

PeriodicMap Nektar::MultiRegions::DisContField2D::m_periodicEdges
protected

A map which identifies pairs of periodic edges.

Definition at line 137 of file DisContField2D.h.

Referenced by Nektar::MultiRegions::ContField2D::ContField2D(), DisContField2D(), FindPeriodicEdges(), IsLeftAdjacentEdge(), SetUpDG(), and v_GetPeriodicEntities().

vector<int> Nektar::MultiRegions::DisContField2D::m_periodicFwdCopy
protected

A vector indicating degress of freedom which need to be copied from forwards to backwards space in case of a periodic boundary condition.

Definition at line 145 of file DisContField2D.h.

Referenced by DisContField2D(), SetUpDG(), and v_GetFwdBwdTracePhys().

PeriodicMap Nektar::MultiRegions::DisContField2D::m_periodicVerts
protected

A map which identifies groups of periodic vertices.

Definition at line 132 of file DisContField2D.h.

Referenced by Nektar::MultiRegions::ContField2D::ContField2D(), DisContField2D(), FindPeriodicEdges(), and v_GetPeriodicEntities().

Array<OneD, Array<OneD, unsigned int> > Nektar::MultiRegions::DisContField2D::m_signEdgeToElmn
protected

Definition at line 121 of file DisContField2D.h.

ExpListSharedPtr Nektar::MultiRegions::DisContField2D::m_trace
protected
AssemblyMapDGSharedPtr Nektar::MultiRegions::DisContField2D::m_traceMap
protected