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

#include <DisContField3D.h>

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

Public Member Functions

 DisContField3D ()
 Default constructor.
 DisContField3D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &graph3D, const std::string &variable, const bool SetUpJustDG=true)
 Constructs a global discontinuous field based on an input mesh with boundary conditions.
 DisContField3D (const DisContField3D &In, const SpatialDomains::MeshGraphSharedPtr &graph3D, const std::string &variable, const bool SetUpJustDG=false)
 DisContField3D (const DisContField3D &In)
 Constructs a global discontinuous field based on another discontinuous field.
virtual ~DisContField3D ()
 Destructor.
GlobalLinSysSharedPtr GetGlobalBndLinSys (const GlobalLinSysKey &mkey)
void EvaluateHDGPostProcessing (Array< OneD, NekDouble > &outarray)
 Evaluate HDG post-processing to increase polynomial order of solution.
- Public Member Functions inherited from Nektar::MultiRegions::ExpList3D
 ExpList3D ()
 Default constructor.
 ExpList3D (const ExpList3D &In)
 Copy constructor.
 ExpList3D (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &TBa, const LibUtilities::BasisKey &TBb, const LibUtilities::BasisKey &TBc, const LibUtilities::BasisKey &HBa, const LibUtilities::BasisKey &HBb, const LibUtilities::BasisKey &HBc, const SpatialDomains::MeshGraphSharedPtr &graph3D, const LibUtilities::PointsType TetNb=LibUtilities::SIZE_PointsType)
 ExpList3D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &graph3D, const std::string &variable="DefaultVar")
 Sets up a list of local expansions based on an input mesh.
 ExpList3D (const SpatialDomains::ExpansionMap &expansions)
 Sets up a list of local expansions based on an expansion vector.
virtual ~ExpList3D ()
 Destructor.
- Public Member Functions inherited from Nektar::MultiRegions::ExpList
 ExpList ()
 The default constructor.
 ExpList (const LibUtilities::SessionReaderSharedPtr &pSession)
 The default constructor.
 ExpList (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 The default constructor.
 ExpList (const ExpList &in, const bool DeclareCoeffPhysArrays=true)
 The copy constructor.
virtual ~ExpList ()
 The default destructor.
int GetNcoeffs (void) const
 Returns the total number of local degrees of freedom $N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_m$.
int GetNcoeffs (const int eid) const
 Returns the total number of local degrees of freedom for element eid.
ExpansionType GetExpType (void)
 Returns the type of the expansion.
void SetExpType (ExpansionType Type)
 Returns the type of the expansion.
int EvalBasisNumModesMax (void) const
 Evaulates the maximum number of modes in the elemental basis order over all elements.
const Array< OneD, int > EvalBasisNumModesMaxPerExp (void) const
 Returns the vector of the number of modes in the elemental basis order over all elements.
int GetTotPoints (void) const
 Returns the total number of quadrature points m_npoints $=Q_{\mathrm{tot}}$.
int GetTotPoints (const int eid) const
 Returns the total number of quadrature points for eid's element $=Q_{\mathrm{tot}}$.
int GetNpoints (void) const
 Returns the total number of quadrature points m_npoints $=Q_{\mathrm{tot}}$.
int Get1DScaledTotPoints (const NekDouble scale) const
 Returns the total number of qudature points scaled by the factor scale on each 1D direction.
void SetWaveSpace (const bool wavespace)
 Sets the wave space to the one of the possible configuration true or false.
void SetModifiedBasis (const bool modbasis)
 Set Modified Basis for the stability analysis.
void SetPhys (int i, NekDouble val)
 Set the i th value of m_phys to value val.
bool GetWaveSpace (void) const
 This function returns the third direction expansion condition, which can be in wave space (coefficient) or not It is stored in the variable m_WaveSpace.
void SetPhys (const Array< OneD, const NekDouble > &inarray)
 Fills the array m_phys.
void SetPhysArray (Array< OneD, NekDouble > &inarray)
 Sets the array m_phys.
void SetPhysState (const bool physState)
 This function manually sets whether the array of physical values $\boldsymbol{u}_l$ (implemented as m_phys) is filled or not.
bool GetPhysState (void) const
 This function indicates whether the array of physical values $\boldsymbol{u}_l$ (implemented as m_phys) is filled or not.
NekDouble PhysIntegral (void)
 This function integrates a function $f(\boldsymbol{x})$ over the domain consisting of all the elements of the expansion.
NekDouble PhysIntegral (const Array< OneD, const NekDouble > &inarray)
 This function integrates a function $f(\boldsymbol{x})$ over the domain consisting of all the elements of the expansion.
void IProductWRTBase_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function calculates the inner product of a function $f(\boldsymbol{x})$ with respect to all {local} expansion modes $\phi_n^e(\boldsymbol{x})$.
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function calculates the inner product of a function $f(\boldsymbol{x})$ with respect to the derivative (in direction.
void FwdTrans_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function elementally evaluates the forward transformation of a function $u(\boldsymbol{x})$ onto the global spectral/hp expansion.
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
void MultiplyByElmtInvMass (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function elementally mulplies the coefficient space of Sin my the elemental inverse of the mass matrix.
void MultiplyByInvMassMatrix (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
void SmoothField (Array< OneD, NekDouble > &field)
 Smooth a field across elements.
void HelmSolve (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const FlagList &flags, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff=StdRegions::NullVarCoeffMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 Solve helmholtz problem.
void LinearAdvectionDiffusionReactionSolve (const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 Solve Advection Diffusion Reaction.
void LinearAdvectionReactionSolve (const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 Solve Advection Diffusion Reaction.
void FwdTrans_BndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void BwdTrans_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function elementally evaluates the backward transformation of the global spectral/hp element expansion.
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
void GetCoords (Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1=NullNekDouble1DArray, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
 This function calculates the coordinates of all the elemental quadrature points $\boldsymbol{x}_i$.
void HomogeneousFwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
void HomogeneousBwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
void DealiasedProd (const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
void GetBCValues (Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
void NormVectorIProductWRTBase (Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
void ApplyGeomInfo ()
 Apply geometry information to each expansion.
void WriteTecplotHeader (std::ostream &outfile, std::string var="")
void WriteTecplotZone (std::ostream &outfile, int expansion=-1)
void WriteTecplotField (std::ostream &outfile, int expansion=-1)
void WriteTecplotConnectivity (std::ostream &outfile, int expansion=-1)
void WriteVtkHeader (std::ostream &outfile)
void WriteVtkFooter (std::ostream &outfile)
void WriteVtkPieceHeader (std::ostream &outfile, int expansion)
void WriteVtkPieceFooter (std::ostream &outfile, int expansion)
void WriteVtkPieceData (std::ostream &outfile, int expansion, std::string var="v")
int GetCoordim (int eid)
 This function returns the dimension of the coordinates of the element eid.
void SetCoeff (int i, NekDouble val)
 Set the i th coefficiient in m_coeffs to value val.
void SetCoeffs (int i, NekDouble val)
 Set the i th coefficiient in m_coeffs to value val.
void SetCoeffsArray (Array< OneD, NekDouble > &inarray)
 Set the m_coeffs array to inarray.
const Array< OneD, const
NekDouble > & 
GetCoeffs () const
 This function returns (a reference to) the array $\boldsymbol{\hat{u}}_l$ (implemented as m_coeffs) containing all local expansion coefficients.
void ImposeDirichletConditions (Array< OneD, NekDouble > &outarray)
 Impose Dirichlet Boundary Conditions onto Array.
void FillBndCondFromField (void)
 Fill Bnd Condition expansion from the values stored in expansion.
void LocalToGlobal (void)
 Put the coefficients into global ordering using m_coeffs.
void GlobalToLocal (void)
 Put the coefficients into local ordering and place in m_coeffs.
NekDouble GetCoeff (int i)
 Get the i th value (coefficient) of m_coeffs.
NekDouble GetCoeffs (int i)
 Get the i th value (coefficient) of m_coeffs.
const Array< OneD, const
NekDouble > & 
GetPhys () const
 This function returns (a reference to) the array $\boldsymbol{u}_l$ (implemented as m_phys) containing the function $u^{\delta}(\boldsymbol{x})$ evaluated at the quadrature points.
NekDouble Linf (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
 This function calculates the $L_\infty$ error of the global spectral/hp element approximation.
NekDouble L2 (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
 This function calculates the $L_2$ error with respect to soln of the global spectral/hp element approximation.
NekDouble H1 (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
 Calculates the $H^1$ error of the global spectral/hp element approximation.
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
Array< OneD, const NekDoubleHomogeneousEnergy (void)
 This function calculates the energy associated with each one of the modesof a 3D homogeneous nD expansion.
void SetHomo1DSpecVanVisc (Array< OneD, NekDouble > visc)
 This function sets the Spectral Vanishing Viscosity in homogeneous1D expansion.
Array< OneD, const unsigned int > GetZIDs (void)
 This function returns a vector containing the wave numbers in z-direction associated with the 3D homogenous expansion. Required if a parellelisation is applied in the Fourier direction.
LibUtilities::TranspositionSharedPtr GetTransposition (void)
 This function returns the transposition class associaed with the homogeneous expansion.
NekDouble GetHomoLen (void)
 This function returns the Width of homogeneous direction associaed with the homogeneous expansion.
Array< OneD, const unsigned int > GetYIDs (void)
 This function returns a vector containing the wave numbers in y-direction associated with the 3D homogenous expansion. Required if a parellelisation is applied in the Fourier direction.
void PhysInterp1DScaled (const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function interpolates the physical space points in inarray to outarray using the same points defined in the expansion but where the number of points are rescaled by 1DScale.
void PhysGalerkinProjection1DScaled (const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function Galerkin projects the physical space points in inarray to outarray where inarray is assumed to be defined in the expansion but where the number of points are rescaled by 1DScale.
int GetExpSize (void)
 This function returns the number of elements in the expansion.
int GetNumElmts (void)
 This function returns the number of elements in the expansion which may be different for a homogeoenous extended expansionp.
const boost::shared_ptr
< LocalRegions::ExpansionVector
GetExp () const
 This function returns the vector of elements in the expansion.
LocalRegions::ExpansionSharedPtrGetExp (int n) const
 This function returns (a shared pointer to) the local elemental expansion of the $n^{\mathrm{th}}$ element.
LocalRegions::ExpansionSharedPtrGetExp (const Array< OneD, const NekDouble > &gloCoord)
 This function returns (a shared pointer to) the local elemental expansion containing the arbitrary point given by gloCoord.
int GetExpIndex (const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0, bool returnNearestElmt=false)
int GetExpIndex (const Array< OneD, const NekDouble > &gloCoords, Array< OneD, NekDouble > &locCoords, NekDouble tol=0.0, bool returnNearestElmt=false)
int GetCoeff_Offset (int n) const
 Get the start offset position for a global list of m_coeffs correspoinding to element n.
int GetPhys_Offset (int n) const
 Get the start offset position for a global list of m_phys correspoinding to element n.
int GetOffset_Elmt_Id (int n) const
 Get the element id associated with the n th consecutive block of data in m_phys and m_coeffs.
Array< OneD, NekDouble > & UpdateCoeffs ()
 This function returns (a reference to) the array $\boldsymbol{\hat{u}}_l$ (implemented as m_coeffs) containing all local expansion coefficients.
Array< OneD, NekDouble > & UpdatePhys ()
 This function returns (a reference to) the array $\boldsymbol{u}_l$ (implemented as m_phys) containing the function $u^{\delta}(\boldsymbol{x})$ evaluated at the quadrature points.
void PhysDeriv (Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 This function discretely evaluates the derivative of a function $f(\boldsymbol{x})$ on the domain consisting of all elements of the expansion.
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
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 DisContField3D &In)
void GenerateBoundaryConditionExpansion (const SpatialDomains::MeshGraphSharedPtr &graph3D, const SpatialDomains::BoundaryConditions &bcs, const std::string &variable)
void FindPeriodicFaces (const SpatialDomains::BoundaryConditions &bcs, const std::string &variable)
 Determine the periodic faces, edges and vertices for the given graph.
bool IsLeftAdjacentFace (const int n, const int e)
virtual void v_GetFwdBwdTracePhys (Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
 This method extracts the "forward" and "backward" trace data from the array field and puts the data into output vectors Fwd and Bwd.
virtual void v_GetFwdBwdTracePhys (const Array< OneD, const NekDouble > &field, Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
virtual void v_ExtractTracePhys (Array< OneD, NekDouble > &outarray)
virtual void v_ExtractTracePhys (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_AddTraceIntegral (const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
 Add trace contributions into elemental coefficient spaces.
virtual void v_AddFwdBwdTraceIntegral (const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
 Add trace contributions into elemental coefficient spaces.
virtual void v_HelmSolve (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const FlagList &flags, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff, const Array< OneD, const NekDouble > &dirForcing)
virtual void v_GeneralMatrixOp (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
 Calculates the result of the multiplication of a global matrix of type specified by mkey with a vector given by inarray.
virtual void v_GetBoundaryToElmtMap (Array< OneD, int > &ElmtID, Array< OneD, int > &FaceID)
 Set up a list of elemeent IDs and edge IDs that link to the boundary conditions.
virtual void v_GetPeriodicEntities (PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces)
virtual ExpListSharedPtrv_GetTrace ()
virtual AssemblyMapDGSharedPtrv_GetTraceMap ()
virtual const Array< OneD,
const
MultiRegions::ExpListSharedPtr > & 
v_GetBndCondExpansions ()
virtual const Array< OneD,
const
SpatialDomains::BoundaryConditionShPtr > & 
v_GetBndConditions ()
virtual
MultiRegions::ExpListSharedPtr
v_UpdateBndCondExpansion (int i)
virtual Array< OneD,
SpatialDomains::BoundaryConditionShPtr > & 
v_UpdateBndConditions ()
virtual void v_EvaluateBoundaryConditions (const NekDouble time=0.0, const std::string varName="", const NekDouble x2_in=NekConstants::kNekUnsetDouble, const NekDouble x3_in=NekConstants::kNekUnsetDouble)
 This function evaluates the boundary conditions at a certain time-level.
virtual map< int,
RobinBCInfoSharedPtr
v_GetRobinBCInfo ()
- Protected Member Functions inherited from Nektar::MultiRegions::ExpList3D
virtual void v_SetUpPhysNormals ()
 Set up the normals on each expansion.
- Protected Member Functions inherited from Nektar::MultiRegions::ExpList
boost::shared_ptr< DNekMatGenGlobalMatrixFull (const GlobalLinSysKey &mkey, const boost::shared_ptr< AssemblyMapCG > &locToGloMap)
const DNekScalBlkMatSharedPtr GenBlockMatrix (const GlobalMatrixKey &gkey)
 This function assembles the block diagonal matrix of local matrices of the type mtype.
const DNekScalBlkMatSharedPtrGetBlockMatrix (const GlobalMatrixKey &gkey)
void MultiplyByBlockMatrix (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
boost::shared_ptr< GlobalMatrixGenGlobalMatrix (const GlobalMatrixKey &mkey, const boost::shared_ptr< AssemblyMapCG > &locToGloMap)
 Generates a global matrix from the given key and map.
void GlobalEigenSystem (const boost::shared_ptr< DNekMat > &Gmat, Array< OneD, NekDouble > &EigValsReal, Array< OneD, NekDouble > &EigValsImag, Array< OneD, NekDouble > &EigVecs=NullNekDouble1DArray)
boost::shared_ptr< GlobalLinSysGenGlobalLinSys (const GlobalLinSysKey &mkey, const boost::shared_ptr< AssemblyMapCG > &locToGloMap)
 This operation constructs the global linear system of type mkey.
boost::shared_ptr< GlobalLinSysGenGlobalBndLinSys (const GlobalLinSysKey &mkey, const AssemblyMapSharedPtr &locToGloMap)
 Generate a GlobalLinSys from information provided by the key "mkey" and the mapping provided in LocToGloBaseMap.
void ReadGlobalOptimizationParameters ()
virtual int v_GetNumElmts (void)
virtual void v_Upwind (const Array< OneD, const Array< OneD, NekDouble > > &Vec, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
virtual void v_Upwind (const Array< OneD, const NekDouble > &Vn, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
virtual const Array< OneD,
const int > & 
v_GetTraceBndMap ()
virtual void v_GetNormals (Array< OneD, Array< OneD, NekDouble > > &normals)
virtual void v_AddTraceIntegral (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
virtual void v_MultiplyByInvMassMatrix (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
virtual void v_LinearAdvectionDiffusionReactionSolve (const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
virtual void v_LinearAdvectionReactionSolve (const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
virtual void v_ImposeDirichletConditions (Array< OneD, NekDouble > &outarray)
virtual void v_FillBndCondFromField ()
virtual void v_LocalToGlobal (void)
virtual void v_GlobalToLocal (void)
virtual void v_BwdTrans (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_ReadGlobalOptimizationParameters ()
virtual std::vector
< LibUtilities::FieldDefinitionsSharedPtr
v_GetFieldDefinitions (void)
virtual void v_GetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
virtual void v_AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
virtual void v_AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, Array< OneD, NekDouble > &coeffs)
virtual void v_ExtractDataToCoeffs (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
 Extract data from raw field data into expansion list.
virtual void v_ExtractCoeffsToCoeffs (const boost::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
virtual void v_WriteTecplotHeader (std::ostream &outfile, std::string var="")
virtual void v_WriteTecplotZone (std::ostream &outfile, int expansion)
virtual void v_WriteTecplotField (std::ostream &outfile, int expansion)
virtual void v_WriteTecplotConnectivity (std::ostream &outfile, int expansion)
virtual void v_WriteVtkPieceHeader (std::ostream &outfile, int expansion)
virtual void v_WriteVtkPieceData (std::ostream &outfile, int expansion, std::string var)
virtual NekDouble v_L2 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray)
virtual Array< OneD, const
NekDouble
v_HomogeneousEnergy (void)
virtual
LibUtilities::TranspositionSharedPtr 
v_GetTransposition (void)
virtual NekDouble v_GetHomoLen (void)
virtual Array< OneD, const
unsigned int > 
v_GetZIDs (void)
virtual Array< OneD, const
unsigned int > 
v_GetYIDs (void)
virtual void v_PhysInterp1DScaled (const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_PhysGalerkinProjection1DScaled (const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void ExtractFileBCs (const std::string &fileName, const std::string &varName, const boost::shared_ptr< ExpList > locExpList)

Protected Attributes

Array< OneD,
MultiRegions::ExpListSharedPtr
m_bndCondExpansions
 An object which contains the discretised boundary conditions.
Array< OneD,
SpatialDomains::BoundaryConditionShPtr
m_bndConditions
 An array which contains the information about the boundary condition on the different boundary regions.
GlobalLinSysMapShPtr m_globalBndMat
ExpListSharedPtr m_trace
AssemblyMapDGSharedPtr m_traceMap
std::set< int > m_boundaryFaces
 A set storing the global IDs of any boundary faces.
PeriodicMap m_periodicFaces
 A map which identifies pairs of periodic faces.
PeriodicMap m_periodicEdges
 A map which identifies groups of periodic edges.
PeriodicMap m_periodicVerts
 A map which identifies groups of periodic vertices.
vector< bool > m_leftAdjacentFaces
vector< int > m_periodicFwdCopy
 A vector indicating degress of freedom which need to be copied from forwards to backwards space in case of a periodic boundary condition.
vector< int > m_periodicBwdCopy

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

Abstraction of a global discontinuous three-dimensional spectral/hp element expansion which approximates the solution of a set of partial differential equations.

Definition at line 51 of file DisContField3D.h.

Constructor & Destructor Documentation

Nektar::MultiRegions::DisContField3D::DisContField3D ( )

Default constructor.

Definition at line 63 of file DisContField3D.cpp.

Nektar::MultiRegions::DisContField3D::DisContField3D ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr graph3D,
const std::string &  variable,
const bool  SetUpJustDG = true 
)

Constructs a global discontinuous field based on an input mesh with boundary conditions.

Definition at line 75 of file DisContField3D.cpp.

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

: ExpList3D (pSession,graph3D,variable),
{
if(variable.compare("DefaultVar") != 0) // do not set up BCs if default variable
{
GenerateBoundaryConditionExpansion(graph3D,bcs,variable);
// Find periodic edges for this variable.
FindPeriodicFaces(bcs, variable);
}
if(SetUpJustDG)
{
}
else
{
// Set element edges to point to Robin BC edges if required.
int i,cnt,f;
Array<OneD, int> ElmtID, FaceID;
GetBoundaryToElmtMap(ElmtID, FaceID);
for(cnt = i = 0; i < m_bndCondExpansions.num_elements(); ++i)
{
locExpList = m_bndCondExpansions[i];
for(f = 0; f < locExpList->GetExpSize(); ++f)
{
= (*m_exp)[ElmtID[cnt+f]]->
as<LocalRegions::Expansion3D>();
= locExpList->GetExp(f)->
as<LocalRegions::Expansion2D>();
exp3d->SetFaceExp(FaceID[cnt+f],exp2d);
exp2d->SetAdjacentElementExp(FaceID[cnt+f],exp3d);
}
cnt += m_bndCondExpansions[i]->GetExpSize();
}
}
}
Nektar::MultiRegions::DisContField3D::DisContField3D ( const DisContField3D In,
const SpatialDomains::MeshGraphSharedPtr graph3D,
const std::string &  variable,
const bool  SetUpJustDG = false 
)

Definition at line 133 of file DisContField3D.cpp.

References Nektar::MultiRegions::ExpList::ApplyGeomInfo(), Nektar::MultiRegions::ExpList::EvaluateBoundaryConditions(), FindPeriodicFaces(), GenerateBoundaryConditionExpansion(), Nektar::MultiRegions::ExpList::GetBoundaryToElmtMap(), m_bndCondExpansions, m_globalBndMat, m_periodicEdges, m_periodicFaces, m_periodicVerts, Nektar::MultiRegions::ExpList::m_session, m_trace, m_traceMap, SameTypeOfBoundaryConditions(), SetUpDG(), and Nektar::MultiRegions::ExpList::SetUpPhysNormals().

: ExpList3D(In),
{
GenerateBoundaryConditionExpansion(graph3D,bcs,variable);
{
// Find periodic edges for this variable.
FindPeriodicFaces(bcs, variable);
if (SetUpJustDG)
{
SetUpDG(variable);
}
else
{
int i,cnt,f;
Array<OneD, int> ElmtID,FaceID;
GetBoundaryToElmtMap(ElmtID,FaceID);
for(cnt = i = 0; i < m_bndCondExpansions.num_elements(); ++i)
{
locExpList = m_bndCondExpansions[i];
for(f = 0; f < locExpList->GetExpSize(); ++f)
{
= (*m_exp)[ElmtID[cnt+f]]->
as<LocalRegions::Expansion3D>();
= locExpList->GetExp(f)->
as<LocalRegions::Expansion2D>();
exp3d->SetFaceExp(FaceID[cnt+f],exp2d);
exp2d->SetAdjacentElementExp(FaceID[cnt+f],exp3d);
}
cnt += m_bndCondExpansions[i]->GetExpSize();
}
}
}
//else if we have the same boundary condition
else
{
if(SetUpJustDG)
{
m_globalBndMat = In.m_globalBndMat;
m_trace = In.m_trace;
m_traceMap = In.m_traceMap;
m_periodicVerts = In.m_periodicVerts;
m_periodicEdges = In.m_periodicEdges;
m_periodicFaces = In.m_periodicFaces;
}
else
{
m_globalBndMat = In.m_globalBndMat;
m_trace = In.m_trace;
m_traceMap = In.m_traceMap;
m_periodicVerts = In.m_periodicVerts;
m_periodicEdges = In.m_periodicEdges;
m_periodicFaces = In.m_periodicFaces;
int i,cnt,f;
Array<OneD, int> ElmtID,FaceID;
GetBoundaryToElmtMap(ElmtID,FaceID);
for(cnt = i = 0; i < m_bndCondExpansions.num_elements(); ++i)
{
locExpList = m_bndCondExpansions[i];
for(f = 0; f < locExpList->GetExpSize(); ++f)
{
= (*m_exp)[ElmtID[cnt+f]]->
as<LocalRegions::Expansion3D>();
= locExpList->GetExp(f)->
as<LocalRegions::Expansion2D>();
exp3d->SetFaceExp(FaceID[cnt+f],exp2d);
exp2d->SetAdjacentElementExp(FaceID[cnt+f],exp3d);
}
cnt += m_bndCondExpansions[i]->GetExpSize();
}
if(m_session->DefinesSolverInfo("PROJECTION"))
{
std::string ProjectStr =
m_session->GetSolverInfo("PROJECTION");
if (ProjectStr == "MixedCGDG" ||
ProjectStr == "Mixed_CG_Discontinuous")
{
SetUpDG(variable);
}
else
{
}
}
else
{
}
}
}
}
Nektar::MultiRegions::DisContField3D::DisContField3D ( const DisContField3D In)

Constructs a global discontinuous field based on another discontinuous field.

Definition at line 257 of file DisContField3D.cpp.

:
ExpList3D(In),
m_bndCondExpansions (In.m_bndCondExpansions),
m_bndConditions (In.m_bndConditions),
m_globalBndMat (In.m_globalBndMat),
m_trace (In.m_trace),
m_traceMap (In.m_traceMap),
m_periodicFaces (In.m_periodicFaces),
m_periodicEdges (In.m_periodicEdges),
m_periodicVerts (In.m_periodicVerts)
{
}
Nektar::MultiRegions::DisContField3D::~DisContField3D ( )
virtual

Destructor.

Definition at line 273 of file DisContField3D.cpp.

{
}

Member Function Documentation

void Nektar::MultiRegions::DisContField3D::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 2236 of file DisContField3D.cpp.

References ASSERTL0, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauMAlpha2Beta0, Nektar::LibUtilities::eHexahedron, Nektar::StdRegions::eInvLaplacianWithUnityMean, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrtho_B, Nektar::LibUtilities::eOrtho_C, Nektar::LibUtilities::ePrism, Nektar::LibUtilities::eTetrahedron, 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,f,ncoeff_face;
Array<OneD, NekDouble> force, out_tmp,qrhs,qrhs1;
Array<OneD, Array< OneD, LocalRegions::ExpansionSharedPtr> >
&elmtToTrace = m_traceMap->GetElmtToTrace();
int eid,nq_elmt, nm_elmt;
int LocBndCoeffs = m_traceMap->GetNumLocalBndCoeffs();
Array<OneD, NekDouble> loc_lambda(LocBndCoeffs), face_lambda;
Array<OneD, NekDouble> tmp_coeffs;
m_traceMap->GlobalToLocalBnd(m_trace->GetCoeffs(),loc_lambda);
face_lambda = loc_lambda;
// Calculate Q using standard DG formulation.
for(i = cnt = 0; i < GetExpSize(); ++i)
{
(*m_exp)[i]->as<LocalRegions::Expansion3D>();
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_points2 = (*m_exp)[eid]->GetBasis(2)->GetNumPoints();
int num_modes0 = (*m_exp)[eid]->GetBasis(0)->GetNumModes();
int num_modes1 = (*m_exp)[eid]->GetBasis(1)->GetNumModes();
int num_modes2 = (*m_exp)[eid]->GetBasis(2)->GetNumModes();
// Probably a better way of setting up lambda than this. Note
// cannot use PutCoeffsInToElmts since lambda space is mapped
// during the solve.
int nFaces = (*m_exp)[eid]->GetNfaces();
Array<OneD, Array<OneD, NekDouble> > faceCoeffs(nFaces);
for(f = 0; f < nFaces; ++f)
{
ncoeff_face = elmtToTrace[eid][f]->GetNcoeffs();
faceCoeffs[f] = Array<OneD, NekDouble>(ncoeff_face);
Vmath::Vcopy(ncoeff_face, face_lambda, 1, faceCoeffs[f], 1);
exp->SetFaceToGeomOrientation(f, faceCoeffs[f]);
face_lambda = face_lambda + ncoeff_face;
}
//creating orthogonal expansion (checking if we have quads or triangles)
LibUtilities::ShapeType shape = (*m_exp)[eid]->DetShapeType();
switch(shape)
{
{
LibUtilities::BasisKey BkeyH1(LibUtilities::eOrtho_A, num_modes0, PkeyH1);
LibUtilities::BasisKey BkeyH2(LibUtilities::eOrtho_A, num_modes1, PkeyH2);
LibUtilities::BasisKey BkeyH3(LibUtilities::eOrtho_A, num_modes2, PkeyH3);
SpatialDomains::HexGeomSharedPtr hGeom = boost::dynamic_pointer_cast<SpatialDomains::HexGeom>((*m_exp)[eid]->GetGeom());
ppExp = MemoryManager<LocalRegions::HexExp>::AllocateSharedPtr(BkeyH1, BkeyH2, BkeyH3, hGeom);
}
break;
{
LibUtilities::BasisKey BkeyT1(LibUtilities::eOrtho_A, num_modes0, PkeyT1);
LibUtilities::BasisKey BkeyT2(LibUtilities::eOrtho_B, num_modes1, PkeyT2);
LibUtilities::BasisKey BkeyT3(LibUtilities::eOrtho_C, num_modes2, PkeyT3);
SpatialDomains::TetGeomSharedPtr tGeom = boost::dynamic_pointer_cast<SpatialDomains::TetGeom>((*m_exp)[eid]->GetGeom());
ppExp = MemoryManager<LocalRegions::TetExp>::AllocateSharedPtr(BkeyT1, BkeyT2, BkeyT3, tGeom);
}
break;
{
LibUtilities::BasisKey BkeyP1(LibUtilities::eOrtho_A, num_modes0, PkeyP1);
LibUtilities::BasisKey BkeyP2(LibUtilities::eOrtho_A, num_modes1, PkeyP2);
LibUtilities::BasisKey BkeyP3(LibUtilities::eOrtho_B, num_modes2, PkeyP3);
SpatialDomains::PrismGeomSharedPtr pGeom = boost::dynamic_pointer_cast<SpatialDomains::PrismGeom>((*m_exp)[eid]->GetGeom());
ppExp = MemoryManager<LocalRegions::PrismExp>::AllocateSharedPtr(BkeyP1, BkeyP2, BkeyP3, pGeom);
}
break;
default:
ASSERTL0(false, "Wrong shape type, HDG postprocessing is not implemented");
};
//DGDeriv
// (d/dx w, q_0)
(*m_exp)[eid]->DGDeriv(
0,tmp_coeffs = m_coeffs + m_coeff_offset[eid],
elmtToTrace[eid], faceCoeffs, out_tmp);
(*m_exp)[eid]->BwdTrans(out_tmp,qrhs);
ppExp->IProductWRTDerivBase(0,qrhs,force);
// + (d/dy w, q_1)
(*m_exp)[eid]->DGDeriv(
1,tmp_coeffs = m_coeffs + m_coeff_offset[eid],
elmtToTrace[eid], faceCoeffs, out_tmp);
(*m_exp)[eid]->BwdTrans(out_tmp,qrhs);
ppExp->IProductWRTDerivBase(1,qrhs,out_tmp);
Vmath::Vadd(nm_elmt,force,1,out_tmp,1,force,1);
// + (d/dz w, q_2)
(*m_exp)[eid]->DGDeriv(
2,tmp_coeffs = m_coeffs + m_coeff_offset[eid],
elmtToTrace[eid], faceCoeffs, out_tmp);
(*m_exp)[eid]->BwdTrans(out_tmp,qrhs);
ppExp->IProductWRTDerivBase(2,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::DisContField3D::FindPeriodicFaces ( const SpatialDomains::BoundaryConditions bcs,
const std::string &  variable 
)
protected

Determine the periodic faces, edges and vertices for the given graph.

Parameters
bcsInformation about the boundary conditions.
variableSpecifies the field.

Definition at line 661 of file DisContField3D.cpp.

References ASSERTL0, ASSERTL1, Nektar::StdRegions::eBackwards, Nektar::StdRegions::eDir1BwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1BwdDir1_Dir2FwdDir2, Nektar::StdRegions::eDir1BwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1BwdDir2_Dir2FwdDir1, Nektar::StdRegions::eDir1FwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1FwdDir1_Dir2FwdDir2, Nektar::StdRegions::eDir1FwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1FwdDir2_Dir2FwdDir1, Nektar::StdRegions::eForwards, Nektar::StdRegions::eNoOrientation, Nektar::SpatialDomains::ePeriodic, Nektar::MultiRegions::ExpList::GetBoundaryCondition(), Nektar::SpatialDomains::BoundaryConditions::GetBoundaryConditions(), Nektar::SpatialDomains::BoundaryConditions::GetBoundaryRegions(), Nektar::SpatialDomains::QuadGeom::GetFaceOrientation(), Nektar::SpatialDomains::TriGeom::GetFaceOrientation(), Nektar::iterator, Nektar::MultiRegions::ExpList::m_graph, m_periodicEdges, m_periodicFaces, m_periodicVerts, Nektar::MultiRegions::ExpList::m_session, Nektar::LibUtilities::ReduceSum, and Vmath::Vsum().

Referenced by DisContField3D().

{
= boost::dynamic_pointer_cast<
SpatialDomains::BoundaryRegionCollection::const_iterator it;
m_session->GetComm()->GetRowComm();
m_session->GetCompositeOrdering();
m_session->GetBndRegionOrdering();
m_graph->GetComposites();
// perComps: Stores a 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).
//
// The four maps allVerts, allCoord, allEdges and allOrient map a
// periodic face to a vector containing the vertex ids of the face;
// their coordinates; the edge ids of the face; and their
// orientation within that face respectively.
//
// Finally the three sets locVerts, locEdges and locFaces store any
// vertices, edges and faces that belong to a periodic composite and
// lie on this process.
map<int,int> perComps;
map<int,vector<int> > allVerts;
map<int,SpatialDomains::PointGeomVector> allCoord;
map<int,vector<int> > allEdges;
map<int,vector<StdRegions::Orientation> > allOrient;
set<int> locVerts;
set<int> locEdges;
set<int> locFaces;
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);
}
for(j = 0; j < (*m_exp)[i]->GetNedges(); ++j)
{
int id = (*m_exp)[i]->GetGeom()->GetEid(j);
locEdges.insert(id);
}
}
// Begin by populating the perComps map. We loop over all periodic
// boundary conditions and determine the composite associated with
// it, then fill out the all* maps.
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;
// Check the region only contains a single composite.
ASSERTL0(it->second->size() == 1,
"Boundary region "+boost::lexical_cast<string>(
region1ID)+" should only contain 1 composite.");
// From this identify composites by looking at the original
// boundary region ordering. Note that in serial the mesh
// partitioner is not run, so this map will be empty and
// therefore needs to be populated by using the corresponding
// boundary region.
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];
}
SpatialDomains::Composite c = it->second->begin()->second;
vector<unsigned int> tmpOrder;
// From the composite, we now construct the allVerts, allEdges
// and allCoord map so that they can be transferred across
// processors. We also populate the locFaces set to store a
// record of all faces local to this process.
for (i = 0; i < c->size(); ++i)
{
boost::dynamic_pointer_cast<
ASSERTL1(faceGeom, "Unable to cast to shared ptr");
// Get geometry ID of this face and store in locFaces.
int faceId = (*c)[i]->GetGlobalID();
locFaces.insert(faceId);
// 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());
}
// Loop over vertices and edges of the face to populate
// allVerts, allEdges and allCoord maps.
vector<int> vertList, edgeList;
vector<StdRegions::Orientation> orientVec;
for (j = 0; j < faceGeom->GetNumVerts(); ++j)
{
vertList .push_back(faceGeom->GetVid (j));
edgeList .push_back(faceGeom->GetEid (j));
coordVec .push_back(faceGeom->GetVertex(j));
orientVec.push_back(faceGeom->GetEorient(j));
}
allVerts [faceId] = vertList;
allEdges [faceId] = edgeList;
allCoord [faceId] = coordVec;
allOrient[faceId] = orientVec;
}
// In serial, record the composite ordering in compOrder for
// later in the routine.
if (vComm->GetSize() == 1)
{
compOrder[it->second->begin()->first] = tmpOrder;
}
// See if we already have either region1 or region2 stored in
// perComps map already and do a sanity check to ensure regions
// are mutually periodic.
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());
}
}
// The next routines process local face lists to exchange vertices,
// edges and faces.
int n = vComm->GetSize();
int p = vComm->GetRank();
int totFaces;
Array<OneD, int> facecounts(n,0);
Array<OneD, int> vertcounts(n,0);
Array<OneD, int> faceoffset(n,0);
Array<OneD, int> vertoffset(n,0);
// First exchange the number of faces on each process.
facecounts[p] = locFaces.size();
vComm->AllReduce(facecounts, LibUtilities::ReduceSum);
// Set up an offset map to allow us to distribute face IDs to all
// processors.
faceoffset[0] = 0;
for (i = 1; i < n; ++i)
{
faceoffset[i] = faceoffset[i-1] + facecounts[i-1];
}
// Calculate total number of faces.
totFaces = Vmath::Vsum(n, facecounts, 1);
// faceIds holds face IDs for each periodic face. faceVerts holds
// the number of vertices in this face.
Array<OneD, int> faceIds (totFaces, 0);
Array<OneD, int> faceVerts(totFaces, 0);
// Process p writes IDs of its faces into position faceoffset[p] of
// faceIds which allows us to perform an AllReduce to distribute
// information amongst processors.
for (i = 0, sIt = locFaces.begin(); sIt != locFaces.end(); ++sIt)
{
faceIds [faceoffset[p] + i ] = *sIt;
faceVerts[faceoffset[p] + i++] = allVerts[*sIt].size();
}
vComm->AllReduce(faceIds, LibUtilities::ReduceSum);
vComm->AllReduce(faceVerts, LibUtilities::ReduceSum);
// procVerts holds number of vertices (and also edges since each
// face is 2D) on each process.
Array<OneD, int> procVerts(n,0);
int nTotVerts;
// Note if there are no periodic faces at all calling Vsum will
// cause a segfault.
if (totFaces > 0)
{
// Calculate number of vertices on each processor.
nTotVerts = Vmath::Vsum(totFaces, faceVerts, 1);
}
else
{
nTotVerts = 0;
}
for (i = 0; i < n; ++i)
{
if (facecounts[i] > 0)
{
procVerts[i] = Vmath::Vsum(
facecounts[i], faceVerts + faceoffset[i], 1);
}
else
{
procVerts[i] = 0;
}
}
// vertoffset is defined in the same manner as edgeoffset
// beforehand.
vertoffset[0] = 0;
for (i = 1; i < n; ++i)
{
vertoffset[i] = vertoffset[i-1] + procVerts[i-1];
}
// At this point we exchange all vertex IDs, edge IDs and vertex
// coordinates for each face. The coordinates are necessary because
// we need to calculate relative face orientations between periodic
// faces to determined edge and vertex connectivity.
Array<OneD, int> vertIds(nTotVerts, 0);
Array<OneD, int> edgeIds(nTotVerts, 0);
Array<OneD, int> edgeOrt(nTotVerts, 0);
Array<OneD, NekDouble> vertX (nTotVerts, 0.0);
Array<OneD, NekDouble> vertY (nTotVerts, 0.0);
Array<OneD, NekDouble> vertZ (nTotVerts, 0.0);
for (cnt = 0, sIt = locFaces.begin();
sIt != locFaces.end(); ++sIt)
{
for (j = 0; j < allVerts[*sIt].size(); ++j)
{
int vertId = allVerts[*sIt][j];
vertIds[vertoffset[p] + cnt ] = vertId;
vertX [vertoffset[p] + cnt ] = (*allCoord[*sIt][j])(0);
vertY [vertoffset[p] + cnt ] = (*allCoord[*sIt][j])(1);
vertZ [vertoffset[p] + cnt ] = (*allCoord[*sIt][j])(2);
edgeIds[vertoffset[p] + cnt ] = allEdges [*sIt][j];
edgeOrt[vertoffset[p] + cnt++] = allOrient[*sIt][j];
}
}
vComm->AllReduce(vertIds, LibUtilities::ReduceSum);
vComm->AllReduce(vertX, LibUtilities::ReduceSum);
vComm->AllReduce(vertY, LibUtilities::ReduceSum);
vComm->AllReduce(vertZ, LibUtilities::ReduceSum);
vComm->AllReduce(edgeIds, LibUtilities::ReduceSum);
vComm->AllReduce(edgeOrt, LibUtilities::ReduceSum);
// Finally now we have all of this information, we construct maps
// which make accessing the information easier. These are
// conceptually the same as all* maps at the beginning of the
// routine, but now hold information for all periodic vertices.
map<int, vector<int> > vertMap;
map<int, vector<int> > edgeMap;
map<int, SpatialDomains::PointGeomVector> coordMap;
// These final two maps are required for determining the relative
// orientation of periodic edges. vCoMap associates vertex IDs with
// their coordinates, and eIdMap maps an edge ID to the two vertices
// which construct it.
map<int, SpatialDomains::PointGeomSharedPtr> vCoMap;
map<int, pair<int, int> > eIdMap;
for (cnt = i = 0; i < totFaces; ++i)
{
vector<int> edges(faceVerts[i]);
vector<int> verts(faceVerts[i]);
SpatialDomains::PointGeomVector coord(faceVerts[i]);
// Keep track of cnt to enable correct edge vertices to be
// inserted into eIdMap.
int tmp = cnt;
for (j = 0; j < faceVerts[i]; ++j, ++cnt)
{
edges[j] = edgeIds[cnt];
verts[j] = vertIds[cnt];
3, verts[j], vertX[cnt], vertY[cnt], vertZ[cnt]);
vCoMap[vertIds[cnt]] = coord[j];
// Try to insert edge into the eIdMap to avoid re-inserting.
pair<map<int, pair<int, int> >::iterator, bool> testIns =
eIdMap.insert(make_pair(
edgeIds[cnt],
make_pair(vertIds[tmp+j],
vertIds[tmp+((j+1) % faceVerts[i])])));
if (testIns.second == false)
{
continue;
}
// If the edge is reversed with respect to the face, then
// swap the edges so that we have the original ordering of
// the edge in the 3D element. This is necessary to properly
// determine edge orientation.
if ((StdRegions::Orientation)edgeOrt[cnt]
{
swap(testIns.first->second.first,
testIns.first->second.second);
}
}
vertMap [faceIds[i]] = verts;
edgeMap [faceIds[i]] = edges;
coordMap[faceIds[i]] = coord;
}
// 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;
// Store temporary map of periodic vertices which will hold all
// periodic vertices on the entire mesh so that doubly periodic
// vertices/edges can be counted properly across partitions. Local
// vertices/edges are copied into m_periodicVerts and
// m_periodicEdges at the end of the function.
PeriodicMap periodicVerts;
PeriodicMap periodicEdges;
// Construct two maps which determine how vertices and edges of
// faces connect given a specific face orientation. The key of the
// map is the number of vertices in the face, used to determine
// difference between tris and quads.
map<int, map<StdRegions::Orientation, vector<int> > > vmap;
map<int, map<StdRegions::Orientation, vector<int> > > emap;
map<StdRegions::Orientation, vector<int> > quadVertMap;
quadVertMap[StdRegions::eDir1FwdDir1_Dir2FwdDir2] += 0,1,2,3;
quadVertMap[StdRegions::eDir1FwdDir1_Dir2BwdDir2] += 3,2,1,0;
quadVertMap[StdRegions::eDir1BwdDir1_Dir2FwdDir2] += 1,0,3,2;
quadVertMap[StdRegions::eDir1BwdDir1_Dir2BwdDir2] += 2,3,0,1;
quadVertMap[StdRegions::eDir1FwdDir2_Dir2FwdDir1] += 0,3,2,1;
quadVertMap[StdRegions::eDir1FwdDir2_Dir2BwdDir1] += 1,2,3,0;
quadVertMap[StdRegions::eDir1BwdDir2_Dir2FwdDir1] += 3,0,1,2;
quadVertMap[StdRegions::eDir1BwdDir2_Dir2BwdDir1] += 2,1,0,3;
map<StdRegions::Orientation, vector<int> > quadEdgeMap;
quadEdgeMap[StdRegions::eDir1FwdDir1_Dir2FwdDir2] += 0,1,2,3;
quadEdgeMap[StdRegions::eDir1FwdDir1_Dir2BwdDir2] += 2,1,0,3;
quadEdgeMap[StdRegions::eDir1BwdDir1_Dir2FwdDir2] += 0,3,2,1;
quadEdgeMap[StdRegions::eDir1BwdDir1_Dir2BwdDir2] += 2,3,0,1;
quadEdgeMap[StdRegions::eDir1FwdDir2_Dir2FwdDir1] += 3,2,1,0;
quadEdgeMap[StdRegions::eDir1FwdDir2_Dir2BwdDir1] += 1,2,3,0;
quadEdgeMap[StdRegions::eDir1BwdDir2_Dir2FwdDir1] += 3,0,1,2;
quadEdgeMap[StdRegions::eDir1BwdDir2_Dir2BwdDir1] += 1,0,3,2;
map<StdRegions::Orientation, vector<int> > triVertMap;
map<StdRegions::Orientation, vector<int> > triEdgeMap;
vmap[3] = triVertMap;
vmap[4] = quadVertMap;
emap[3] = triEdgeMap;
emap[4] = quadEdgeMap;
map<int,int> allCompPairs;
// Finally we have enough information to populate the periodic
// vertex, edge and face maps. Begin by looping over all pairs of
// periodic composites to determine pairs of periodic faces.
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],
"Neither composite not found on this process!");
// Loop over composite ordering to construct list of all
// periodic faces, 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!");
// Look up composite ordering to determine pairs.
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.");
// Sanity check that the faces are mutually periodic.
if (compPairs.count(eId2) != 0)
{
ASSERTL0(compPairs[eId2] == eId1, "Pairing incorrect");
}
compPairs[eId1] = eId2;
}
// Now that we have all pairs of periodic faces, loop over the
// ones local to this process and populate face/edge/vertex
// maps.
for (pIt = compPairs.begin(); pIt != compPairs.end(); ++pIt)
{
int ids [2] = {pIt->first, pIt->second};
bool local[2] = {locFaces.count(pIt->first) > 0,
locFaces.count(pIt->second) > 0};
ASSERTL0(coordMap.count(ids[0]) > 0 &&
coordMap.count(ids[1]) > 0,
"Unable to find face in coordinate map");
allCompPairs[pIt->first ] = pIt->second;
allCompPairs[pIt->second] = pIt->first;
// Loop up coordinates of the faces, check they have the
// same number of vertices.
= { coordMap[ids[0]], coordMap[ids[1]] };
ASSERTL0(tmpVec[0].size() == tmpVec[1].size(),
"Two periodic faces have different number "
"of vertices!");
// o will store relative orientation of faces. Note that in
// some transpose cases (Dir1FwdDir2_Dir2BwdDir1 and
// Dir1BwdDir1_Dir2FwdDir1) it seems orientation will be
// different going from face1->face2 instead of face2->face1
// (check this).
// Record periodic faces.
for (i = 0; i < 2; ++i)
{
if (!local[i])
{
continue;
}
// Reference to the other face.
int other = (i+1) % 2;
// Calculate relative face orientation.
if (tmpVec[0].size() == 3)
{
tmpVec[i], tmpVec[other]);
}
else
{
tmpVec[i], tmpVec[other]);
}
// Record face ID, orientation and whether other face is
// local.
PeriodicEntity ent(ids [other], o,
local[other]);
m_periodicFaces[ids[i]].push_back(ent);
}
int nFaceVerts = vertMap[ids[0]].size();
// Determine periodic vertices.
for (i = 0; i < 2; ++i)
{
int other = (i+1) % 2;
// Calculate relative face orientation.
if (tmpVec[0].size() == 3)
{
tmpVec[i], tmpVec[other]);
}
else
{
tmpVec[i], tmpVec[other]);
}
if (nFaceVerts == 3)
{
"Unsupported face orientation for face "+
boost::lexical_cast<string>(ids[i]));
}
// Look up vertices for this face.
vector<int> per1 = vertMap[ids[i]];
vector<int> per2 = vertMap[ids[other]];
// tmpMap will hold the pairs of vertices which are
// periodic.
map<int, pair<int, bool> > tmpMap;
map<int, pair<int, bool> >::iterator mIt;
// Use vmap to determine which vertices connect given
// the orientation o.
for (j = 0; j < nFaceVerts; ++j)
{
int v = vmap[nFaceVerts][o][j];
tmpMap[per1[j]] = make_pair(
per2[v], locVerts.count(per2[v]) > 0);
}
// Now loop over tmpMap to associate periodic vertices.
for (mIt = tmpMap.begin(); mIt != tmpMap.end(); ++mIt)
{
PeriodicEntity ent2(mIt->second.first,
mIt->second.second);
// See if this vertex has been recorded already.
PeriodicMap::iterator perIt = periodicVerts.find(
mIt->first);
if (perIt == periodicVerts.end())
{
// Vertex is new - just record this entity as
// usual.
periodicVerts[mIt->first].push_back(ent2);
perIt = periodicVerts.find(mIt->first);
}
else
{
// Vertex is known - loop over the vertices
// inside the record and potentially add vertex
// mIt->second to the list.
for (k = 0; k < perIt->second.size(); ++k)
{
if (perIt->second[k].id == mIt->second.first)
{
break;
}
}
if (k == perIt->second.size())
{
perIt->second.push_back(ent2);
}
}
}
}
// Determine periodic edges. Logic is the same as above, and
// perhaps should be condensed to avoid replication.
for (i = 0; i < 2; ++i)
{
int other = (i+1) % 2;
if (tmpVec[0].size() == 3)
{
tmpVec[i], tmpVec[other]);
}
else
{
tmpVec[i], tmpVec[other]);
}
vector<int> per1 = edgeMap[ids[i]];
vector<int> per2 = edgeMap[ids[other]];
map<int, pair<int, bool> > tmpMap;
map<int, pair<int, bool> >::iterator mIt;
for (j = 0; j < nFaceVerts; ++j)
{
int e = emap[nFaceVerts][o][j];
tmpMap[per1[j]] = make_pair(
per2[e], locEdges.count(per2[e]) > 0);
}
for (mIt = tmpMap.begin(); mIt != tmpMap.end(); ++mIt)
{
// Note we assume orientation of edges is forwards -
// this may be reversed later.
PeriodicEntity ent2(mIt->second.first,
mIt->second.second);
PeriodicMap::iterator perIt = periodicEdges.find(
mIt->first);
if (perIt == periodicEdges.end())
{
periodicEdges[mIt->first].push_back(ent2);
perIt = periodicEdges.find(mIt->first);
}
else
{
for (k = 0; k < perIt->second.size(); ++k)
{
if (perIt->second[k].id == mIt->second.first)
{
break;
}
}
if (k == perIt->second.size())
{
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 < totFaces; ++i)
{
int faceId = faceIds[i];
ASSERTL0(allCompPairs.count(faceId) > 0,
"Unable to find matching periodic face.");
int perFaceId = allCompPairs[faceId];
for (j = 0; j < faceVerts[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. perFaceId is
// the face ID which is periodic with faceId. The logic is
// much the same as the loop above.
= { coordMap[faceId], coordMap[perFaceId] };
int nFaceVerts = tmpVec[0].size();
StdRegions::Orientation o = nFaceVerts == 3 ?
tmpVec[0], tmpVec[1]) :
SpatialDomains::QuadGeom::GetFaceOrientation(
tmpVec[0], tmpVec[1]);
// Use vmap to determine which vertex of the other face
// should be periodic with this one.
int perVertexId = vertMap[perFaceId][vmap[nFaceVerts][o][j]];
PeriodicEntity ent(perVertexId,
locVerts.count(perVertexId) > 0);
periodicVerts[vId].push_back(ent);
}
int eId = edgeIds[cnt];
perId = periodicEdges.find(eId);
if (perId == periodicEdges.end())
{
// This edge is not included in the map. Figure
// out which edge it is supposed to be periodic
// with. perFaceId is the face ID which is
// periodic with faceId. The logic is much the
// same as the loop above.
= { coordMap[faceId], coordMap[perFaceId] };
int nFaceEdges = tmpVec[0].size();
StdRegions::Orientation o = nFaceEdges == 3 ?
tmpVec[0], tmpVec[1]) :
SpatialDomains::QuadGeom::GetFaceOrientation(
tmpVec[0], tmpVec[1]);
// Use emap to determine which edge of the other
// face should be periodic with this one.
int perEdgeId = edgeMap[perFaceId][emap[nFaceEdges][o][j]];
PeriodicEntity ent(perEdgeId,
locEdges.count(perEdgeId) > 0);
periodicEdges[eId].push_back(ent);
}
}
}
// Finally, we must loop over the periodicVerts and periodicEdges
// map to complete connectivity information.
PeriodicMap::iterator perIt, perIt2;
for (perIt = periodicVerts.begin();
perIt != periodicVerts.end(); ++perIt)
{
// For each vertex that is periodic with this one...
for (i = 0; i < perIt->second.size(); ++i)
{
// Find the vertex in the periodicVerts map...
perIt2 = periodicVerts.find(perIt->second[i].id);
ASSERTL0(perIt2 != periodicVerts.end(),
"Couldn't find periodic vertex.");
// Now search through this vertex's list and make sure that
// we have a record of any vertices which aren't in the
// original list.
for (j = 0; j < perIt2->second.size(); ++j)
{
if (perIt2->second[j].id == perIt->first)
{
continue;
}
for (k = 0; k < perIt->second.size(); ++k)
{
if (perIt2->second[j].id == perIt->second[k].id)
{
break;
}
}
if (k == perIt->second.size())
{
perIt->second.push_back(perIt2->second[j]);
}
}
}
}
for (perIt = periodicEdges.begin();
perIt != periodicEdges.end(); ++perIt)
{
for (i = 0; i < perIt->second.size(); ++i)
{
perIt2 = periodicEdges.find(perIt->second[i].id);
ASSERTL0(perIt2 != periodicEdges.end(),
"Couldn't find periodic edge.");
for (j = 0; j < perIt2->second.size(); ++j)
{
if (perIt2->second[j].id == perIt->first)
{
continue;
}
for (k = 0; k < perIt->second.size(); ++k)
{
if (perIt2->second[j].id == perIt->second[k].id)
{
break;
}
}
if (k == perIt->second.size())
{
perIt->second.push_back(perIt2->second[j]);
}
}
}
}
// Loop over periodic edges to determine relative edge orientations.
for (perIt = periodicEdges.begin();
perIt != periodicEdges.end(); perIt++)
{
// Find edge coordinates
map<int, pair<int, int> >::iterator eIt
= eIdMap.find(perIt->first);
*vCoMap[eIt->second.first],
*vCoMap[eIt->second.second]
};
// Loop over each edge, and construct a vector that takes us
// from one vertex to another. Use this to figure out which
// vertex maps to which.
for (i = 0; i < perIt->second.size(); ++i)
{
eIt = eIdMap.find(perIt->second[i].id);
*vCoMap[eIt->second.first],
*vCoMap[eIt->second.second]
};
NekDouble cx = 0.5*(w[0](0)-v[0](0)+w[1](0)-v[1](0));
NekDouble cy = 0.5*(w[0](1)-v[0](1)+w[1](1)-v[1](1));
NekDouble cz = 0.5*(w[0](2)-v[0](2)+w[1](2)-v[1](2));
int vMap[2] = {-1,-1};
for (j = 0; j < 2; ++j)
{
NekDouble x = v[j](0);
NekDouble y = v[j](1);
NekDouble z = v[j](2);
for (k = 0; k < 2; ++k)
{
NekDouble x1 = w[k](0)-cx;
NekDouble y1 = w[k](1)-cy;
NekDouble z1 = w[k](2)-cz;
if (sqrt((x1-x)*(x1-x)+(y1-y)*(y1-y)+(z1-z)*(z1-z))
< 1e-8)
{
vMap[k] = j;
break;
}
}
}
// Sanity check the map.
ASSERTL0(vMap[0] >= 0 && vMap[1] >= 0,
"Unable to align periodic vertices.");
ASSERTL0((vMap[0] == 0 || vMap[0] == 1) &&
(vMap[1] == 0 || vMap[1] == 1) &&
(vMap[0] != vMap[1]),
"Unable to align periodic vertices.");
// If 0 -> 0 then edges are aligned already; otherwise
// reverse the orientation.
if (vMap[0] != 0)
{
perIt->second[i].orient = StdRegions::eBackwards;
}
}
}
// Do one final loop over periodic vertices/edges to remove
// non-local vertices/edges from map.
for (perIt = periodicVerts.begin();
perIt != periodicVerts.end(); ++perIt)
{
if (locVerts.count(perIt->first) > 0)
{
m_periodicVerts.insert(*perIt);
}
}
for (perIt = periodicEdges.begin();
perIt != periodicEdges.end(); ++perIt)
{
if (locEdges.count(perIt->first) > 0)
{
m_periodicEdges.insert(*perIt);
}
}
}
void Nektar::MultiRegions::DisContField3D::GenerateBoundaryConditionExpansion ( const SpatialDomains::MeshGraphSharedPtr graph3D,
const SpatialDomains::BoundaryConditions bcs,
const std::string &  variable 
)
protected

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

Parameters
graph3DA 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 596 of file DisContField3D.cpp.

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

Referenced by DisContField3D().

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

Definition at line 277 of file DisContField3D.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::DisContField3D::IsLeftAdjacentFace ( const int  n,
const int  e 
)
protected

Definition at line 1622 of file DisContField3D.cpp.

References ASSERTL2, Nektar::iterator, m_boundaryFaces, Nektar::MultiRegions::ExpList::m_exp, m_periodicFaces, m_trace, and m_traceMap.

Referenced by SetUpDG(), and v_AddFwdBwdTraceIntegral().

{
m_traceMap->GetElmtToTrace()[n][e]->
as<LocalRegions::Expansion2D>();
int offset = m_trace->GetPhys_Offset(traceEl->GetElmtId());
bool fwd = true;
if (traceEl->GetLeftAdjacentElementFace () == -1 ||
traceEl->GetRightAdjacentElementFace() == -1)
{
// Boundary edge (1 connected element). Do nothing in
// serial.
it = m_boundaryFaces.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_boundaryFaces.end())
{
int traceGeomId = traceEl->GetGeom2D()->GetGlobalID();
traceGeomId);
if (pIt != m_periodicFaces.end() && !pIt->second[0].isLocal)
{
fwd = traceGeomId == min(traceGeomId,pIt->second[0].id);
}
else
{
fwd = m_traceMap->
GetTraceToUniversalMapUnique(offset) >= 0;
}
}
}
else if (traceEl->GetLeftAdjacentElementFace () != -1 &&
traceEl->GetRightAdjacentElementFace() != -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;
}
bool Nektar::MultiRegions::DisContField3D::SameTypeOfBoundaryConditions ( const DisContField3D In)
protected

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

Parameters
InContField3D to compare with.
Returns
true if boundary conditions match.

Definition at line 554 of file DisContField3D.cpp.

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

Referenced by Nektar::MultiRegions::ContField3D::ContField3D(), and DisContField3D().

{
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->AllReduce(vSame, LibUtilities::ReduceMin);
return (vSame == 1);
}
void Nektar::MultiRegions::DisContField3D::SetUpDG ( const std::string  variable = "DefaultVar")
protected

Set up all DG member variables and maps.

Definition at line 306 of file DisContField3D.cpp.

References Nektar::StdRegions::StdExpansion::as(), ASSERTL1, Nektar::StdRegions::eDir1BwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1BwdDir1_Dir2FwdDir2, Nektar::StdRegions::eDir1BwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1BwdDir2_Dir2FwdDir1, Nektar::StdRegions::eDir1FwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1FwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1FwdDir2_Dir2FwdDir1, Nektar::SpatialDomains::ePeriodic, Nektar::LocalRegions::Expansion3D::GetGeom3D(), Nektar::MultiRegions::_PeriodicEntity::id, IsLeftAdjacentFace(), Nektar::MultiRegions::_PeriodicEntity::isLocal, Nektar::iterator, m_bndCondExpansions, m_bndConditions, m_boundaryFaces, Nektar::MultiRegions::ExpList::m_exp, m_globalBndMat, Nektar::MultiRegions::ExpList::m_graph, m_leftAdjacentFaces, m_periodicBwdCopy, m_periodicFaces, m_periodicFwdCopy, Nektar::MultiRegions::ExpList::m_session, m_trace, m_traceMap, Nektar::MultiRegions::NullExpListSharedPtr, Nektar::MultiRegions::_PeriodicEntity::orient, Nektar::LocalRegions::Expansion2D::SetAdjacentElementExp(), and Nektar::MultiRegions::ExpList::SetUpPhysNormals().

Referenced by DisContField3D(), and v_GetTrace().

{
{
return;
}
boost::dynamic_pointer_cast<SpatialDomains::MeshGraph3D>(
// Set up matrix map
// Set up Trace space
bool UseGenSegExp = true;
*m_exp,graph3D, m_periodicFaces, UseGenSegExp);
m_trace = trace;
m_session,graph3D,trace,*this,m_bndCondExpansions,
Array<OneD, Array<OneD, LocalRegions::ExpansionSharedPtr> >
&elmtToTrace = m_traceMap->GetElmtToTrace();
// Scatter trace segments to 3D 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]->GetNfaces(); ++j)
{
(*m_exp)[i]->as<LocalRegions::Expansion3D>();
elmtToTrace[i][j]->as<LocalRegions::Expansion2D>();
exp3d->SetFaceExp (j, exp2d);
exp2d->SetAdjacentElementExp(j, exp3d);
}
}
// Set up physical normals
// Set up information for parallel jobs.
for (int i = 0; i < m_trace->GetExpSize(); ++i)
{
int offset = m_trace->GetPhys_Offset(i);
int traceGeomId = traceEl->GetGeom2D()->GetGlobalID();
traceGeomId);
if (pIt != m_periodicFaces.end() && !pIt->second[0].isLocal)
{
if (traceGeomId != min(pIt->second[0].id, traceGeomId))
{
traceEl->GetLeftAdjacentElementExp()->NegateFaceNormal(
traceEl->GetLeftAdjacentElementFace());
}
}
else if (m_traceMap->GetTraceToUniversalMapUnique(offset) < 0)
{
traceEl->GetLeftAdjacentElementExp()->NegateFaceNormal(
traceEl->GetLeftAdjacentElementFace());
}
}
int cnt, n, e;
// Identify boundary faces
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> > perFaceToExpMap;
boost::unordered_map<int,pair<int,int> >::iterator it2;
cnt = 0;
for (int n = 0; n < m_exp->size(); ++n)
{
exp3d = (*m_exp)[n]->as<LocalRegions::Expansion3D>();
for (int e = 0; e < exp3d->GetNfaces(); ++e, ++cnt)
{
exp3d->GetGeom3D()->GetFid(e));
if (it != m_periodicFaces.end())
{
perFaceToExpMap[it->first] = make_pair(n, e);
}
}
}
// Set up left-adjacent face list.
m_leftAdjacentFaces.resize(cnt);
cnt = 0;
for (int i = 0; i < m_exp->size(); ++i)
{
for (int j = 0; j < (*m_exp)[i]->GetNfaces(); ++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)
{
exp3d = (*m_exp)[n]->as<LocalRegions::Expansion3D>();
for (int e = 0; e < exp3d->GetNfaces(); ++e, ++cnt)
{
int faceGeomId = exp3d->GetGeom3D()->GetFid(e);
int offset = m_trace->GetPhys_Offset(
elmtToTrace[n][e]->GetElmtId());
// Check to see if this face is periodic.
PeriodicMap::iterator it = m_periodicFaces.find(faceGeomId);
if (it != m_periodicFaces.end())
{
const PeriodicEntity &ent = it->second[0];
it2 = perFaceToExpMap.find(ent.id);
if (it2 == perFaceToExpMap.end())
{
if (m_session->GetComm()->GetSize() > 1 &&
!ent.isLocal)
{
continue;
}
else
{
ASSERTL1(false, "Periodic edge not found!");
}
}
"Periodic face in non-forward space?");
int offset2 = m_trace->GetPhys_Offset(
elmtToTrace[it2->second.first][it2->second.second]->
GetElmtId());
// Calculate relative orientations between faces to
// calculate copying map.
int nquad1 = elmtToTrace[n][e]->GetNumPoints(0);
int nquad2 = elmtToTrace[n][e]->GetNumPoints(1);
vector<int> tmpBwd(nquad1*nquad2);
vector<int> tmpFwd(nquad1*nquad2);
{
for (int i = 0; i < nquad2; ++i)
{
for (int j = 0; j < nquad1; ++j)
{
tmpBwd[i*nquad1+j] = offset2 + i*nquad1+j;
tmpFwd[i*nquad1+j] = offset + j*nquad2+i;
}
}
}
else
{
for (int i = 0; i < nquad2; ++i)
{
for (int j = 0; j < nquad1; ++j)
{
tmpBwd[i*nquad1+j] = offset2 + i*nquad1+j;
tmpFwd[i*nquad1+j] = offset + i*nquad1+j;
}
}
}
{
// Reverse x direction
for (int i = 0; i < nquad2; ++i)
{
for (int j = 0; j < nquad1/2; ++j)
{
swap(tmpFwd[i*nquad1 + j],
tmpFwd[i*nquad1 + nquad1-j-1]);
}
}
}
{
// Reverse y direction
for (int j = 0; j < nquad1; ++j)
{
for (int i = 0; i < nquad2/2; ++i)
{
swap(tmpFwd[i*nquad1 + j],
tmpFwd[(nquad2-i-1)*nquad1 + j]);
}
}
}
for (int i = 0; i < nquad1*nquad2; ++i)
{
m_periodicFwdCopy.push_back(tmpFwd[i]);
m_periodicBwdCopy.push_back(tmpBwd[i]);
}
}
}
}
}
void Nektar::MultiRegions::DisContField3D::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{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. The value of q is determined from the routine IsLeftAdjacentFace() which if true we use Fwd else we use Bwd

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

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1919 of file DisContField3D.cpp.

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

{
int e,n,offset, t_offset;
Array<OneD, NekDouble> e_outarray;
Array<OneD, Array<OneD, LocalRegions::ExpansionSharedPtr> >
&elmtToTrace = m_traceMap->GetElmtToTrace();
for(n = 0; n < GetExpSize(); ++n)
{
offset = GetCoeff_Offset(n);
for(e = 0; e < (*m_exp)[n]->GetNfaces(); ++e)
{
t_offset = m_trace->GetPhys_Offset(elmtToTrace[n][e]->GetElmtId());
// Evaluate upwind flux less local edge
{
(*m_exp)[n]->AddFaceNormBoundaryInt(
e, elmtToTrace[n][e], Fwd + t_offset,
e_outarray = outarray+offset);
}
else
{
(*m_exp)[n]->AddFaceNormBoundaryInt(
e, elmtToTrace[n][e], Bwd + t_offset,
e_outarray = outarray+offset);
}
}
}
}
void Nektar::MultiRegions::DisContField3D::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
Expansion3D::AddFaceNormBoundaryInt
Parameters
FnThe trace quantities.
outarrayResulting 3D coefficient space.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1872 of file DisContField3D.cpp.

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

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

This function evaluates the boundary conditions at a certain time-level.

Based on the boundary condition $g(\boldsymbol{x},t)$ evaluated at a given time-level t, this function transforms the boundary conditions onto the coefficients of the (one-dimensional) boundary expansion. Depending on the type of boundary conditions, these expansion coefficients are calculated in different ways:

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

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 2414 of file DisContField3D.cpp.

References ASSERTL0, Nektar::SpatialDomains::eDirichlet, Nektar::SpatialDomains::eNeumann, Nektar::SpatialDomains::eRobin, Nektar::SpatialDomains::eTimeDependent, Nektar::LibUtilities::Equation::Evaluate(), Nektar::MultiRegions::ExpList::ExtractFileBCs(), 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);
locExpList->GetCoords(x0, x1, x2);
if (m_bndConditions[i]->GetBoundaryConditionType()
{
string filebcs = boost::static_pointer_cast<
m_bndConditions[i])->m_filename;
if (filebcs != "")
{
ExtractFileBCs(filebcs, varName, locExpList);
}
else
{
LibUtilities::Equation condition = boost::static_pointer_cast<SpatialDomains::
m_bndConditions[i])->m_dirichletCondition;
condition.Evaluate(x0, x1, x2, time,
locExpList->UpdatePhys());
locExpList->FwdTrans_BndConstrained(
locExpList->GetPhys(),
locExpList->UpdateCoeffs());
}
}
else if (m_bndConditions[i]->GetBoundaryConditionType()
{
LibUtilities::Equation condition = boost::
static_pointer_cast<SpatialDomains::
m_bndConditions[i])->m_neumannCondition;
condition.Evaluate(x0, x1, x2, time,
locExpList->UpdatePhys());
locExpList->IProductWRTBase(locExpList->GetPhys(),
locExpList->UpdateCoeffs());
}
else if (m_bndConditions[i]->GetBoundaryConditionType()
{
LibUtilities::Equation condition = boost::
static_pointer_cast<SpatialDomains::
m_bndConditions[i])->m_robinFunction;
LibUtilities::Equation coeff = boost::
static_pointer_cast<SpatialDomains::
m_bndConditions[i])->m_robinPrimitiveCoeff;
condition.Evaluate(x0, x1, x2, time,
locExpList->UpdatePhys());
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::DisContField3D::v_ExtractTracePhys ( Array< OneD, NekDouble > &  outarray)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1816 of file DisContField3D.cpp.

References ASSERTL1, Nektar::MultiRegions::ExpList::m_phys, and Nektar::MultiRegions::ExpList::m_physState.

{
"Field is not in physical space.");
}
void Nektar::MultiRegions::DisContField3D::v_ExtractTracePhys ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1825 of file DisContField3D.cpp.

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

{
// 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, LocalRegions::ExpansionSharedPtr> >
&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]->GetNfaces(); ++e)
{
offset = m_trace->GetPhys_Offset(elmtToTrace[n][e]->GetElmtId());
(*m_exp)[n]->GetFacePhysVals(e, elmtToTrace[n][e],
inarray + phys_offset,
e_tmp = outarray + offset);
}
}
}
void Nektar::MultiRegions::DisContField3D::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 2151 of file DisContField3D.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::DisContField3D::v_GetBndCondExpansions ( void  )
inlineprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 215 of file DisContField3D.h.

References m_bndCondExpansions.

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

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 222 of file DisContField3D.h.

References m_bndConditions.

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

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

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1957 of file DisContField3D.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::MeshGraph3D>(
// Populate global ID map (takes global geometry ID to local
// expansion list ID).
for (i = 0; i < GetExpSize(); ++i)
{
exp3d = (*m_exp)[i]->as<LocalRegions::Expansion3D>();
globalIdMap[exp3d->GetGeom3D()->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(FaceID.num_elements() != nbcs)
{
FaceID = Array<OneD, int>(nbcs);
}
for(cnt = n = 0; n < m_bndCondExpansions.num_elements(); ++n)
{
for(i = 0; i < m_bndCondExpansions[n]->GetExpSize(); ++i, ++cnt)
{
exp2d = m_bndCondExpansions[n]->GetExp(i)->
as<LocalRegions::Expansion2D>();
// Use face to element map from MeshGraph3D.
graph3D->GetElementsFromFace(exp2d->GetGeom2D());
ElmtID[cnt] = globalIdMap[(*tmp)[0]->
m_Element->GetGlobalID()];
FaceID[cnt] = (*tmp)[0]->m_FaceIndx;
}
}
}
void Nektar::MultiRegions::DisContField3D::v_GetFwdBwdTracePhys ( 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 face of each element and its corresponding face in the trace space using the mapping m_traceMap. The element can either be left-adjacent or right-adjacent to this trace face (see Expansion2D::GetLeftAdjacentElementExp). Boundary faces are always left-adjacent since left-adjacency is populated first.

If the element is left-adjacent we extract the face 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 3D data from which we wish to extract the backward and forward orientated trace/face arrays.
Returns
Updates a NekDouble array Fwd and Bwd

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1699 of file DisContField3D.cpp.

References Nektar::MultiRegions::ExpList::m_phys.

{
}
void Nektar::MultiRegions::DisContField3D::v_GetFwdBwdTracePhys ( const Array< OneD, const NekDouble > &  field,
Array< OneD, NekDouble > &  Fwd,
Array< OneD, NekDouble > &  Bwd 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1706 of file DisContField3D.cpp.

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

{
// Loop over elements and collect forward and backward expansions.
int nexp = GetExpSize();
int cnt, n, e, npts, offset, phys_offset;
Array<OneD,NekDouble> e_tmp;
Array<OneD, Array<OneD, LocalRegions::ExpansionSharedPtr> >
&elmtToTrace = m_traceMap->GetElmtToTrace();
boost::unordered_map<int,pair<int,int> >::iterator it3;
// Zero vectors.
Vmath::Zero(Fwd.num_elements(), Fwd, 1);
Vmath::Zero(Bwd.num_elements(), Bwd, 1);
bool fwd;
for(cnt = n = 0; n < nexp; ++n)
{
exp3d = (*m_exp)[n]->as<LocalRegions::Expansion3D>();
phys_offset = GetPhys_Offset(n);
for(e = 0; e < exp3d->GetNfaces(); ++e, ++cnt)
{
offset = m_trace->GetPhys_Offset(
elmtToTrace[n][e]->GetElmtId());
fwd = m_leftAdjacentFaces[cnt];
if (fwd)
{
exp3d->GetFacePhysVals(e, elmtToTrace[n][e],
field + phys_offset,
e_tmp = Fwd + offset);
}
else
{
exp3d->GetFacePhysVals(e, elmtToTrace[n][e],
field + phys_offset,
e_tmp = Bwd + offset);
}
}
}
// fill boundary conditions into missing elements
int id1,id2 = 0;
cnt = 0;
for(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)->GetTotPoints();
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)->GetTotPoints();
id1 = m_bndCondExpansions[n]->GetPhys_Offset(e);
id2 = m_trace->GetPhys_Offset(
m_traceMap->GetBndCondTraceToGlobalTraceMap(cnt+e));
ASSERTL0((m_bndCondExpansions[n]->GetPhys())[id1] == 0.0,
"method not set up for non-zero Neumann "
"boundary condition");
Vmath::Vcopy(npts,&Fwd[id2],1,&Bwd[id2],1);
}
cnt += e;
}
else
{
ASSERTL0(false, "Method only set up for Dirichlet, Neumann "
"and Robin conditions.");
}
}
// 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);
}
virtual void Nektar::MultiRegions::DisContField3D::v_GetPeriodicEntities ( PeriodicMap periodicVerts,
PeriodicMap periodicEdges,
PeriodicMap periodicFaces 
)
inlineprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 189 of file DisContField3D.h.

References m_periodicEdges, m_periodicFaces, and m_periodicVerts.

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

Search through the edge expansions and identify which ones have Robin/Mixed type boundary conditions. If find a Robin boundary 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
std map containing the robin boundary condition info using a key of the element id

There is a next member to allow for more than one Robin boundary condition per element

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 2180 of file DisContField3D.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,FaceID;
GetBoundaryToElmtMap(ElmtID,FaceID);
for(cnt = i = 0; i < m_bndCondExpansions.num_elements(); ++i)
{
if(m_bndConditions[i]->GetBoundaryConditionType() == SpatialDomains::eRobin)
{
int e,elmtid;
Array<OneD, NekDouble> Array_tmp;
locExpList = m_bndCondExpansions[i];
for(e = 0; e < locExpList->GetExpSize(); ++e)
{
RobinBCInfoSharedPtr rInfo = MemoryManager<RobinBCInfo>::AllocateSharedPtr(FaceID[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::DisContField3D::v_GetTrace ( )
inlineprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 199 of file DisContField3D.h.

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

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

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 209 of file DisContField3D.h.

References m_traceMap.

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

Solving Helmholtz Equation in 3D

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 2017 of file DisContField3D.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);
const DNekScalBlkMatSharedPtr &HDGLamToU = GetBlockMatrix(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() == SpatialDomains::eDirichlet)
{
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);
const DNekScalBlkMatSharedPtr& InvHDGHelm = GetBlockMatrix(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::DisContField3D::v_UpdateBndCondExpansion ( int  i)
inlineprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 228 of file DisContField3D.h.

References m_bndCondExpansions.

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

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 234 of file DisContField3D.h.

References m_bndConditions.

{
}

Member Data Documentation

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

An object which contains the discretised boundary conditions.

It is an array of size equal to the number of boundary regions and consists of entries of the type MultiRegions::ExpList2D. Every entry corresponds to the two-dimensional spectral/hp expansion on a single boundary region. The values of the boundary conditions are stored as the coefficients of the two-dimensional expansion.

Definition at line 92 of file DisContField3D.h.

Referenced by Nektar::MultiRegions::ContField3D::ContField3D(), DisContField3D(), GenerateBoundaryConditionExpansion(), Nektar::MultiRegions::ContField3D::GenerateDirBndCondForcing(), Nektar::MultiRegions::ContField3D::GetBndCondExp(), Nektar::MultiRegions::ContField3D::GetBndCondExpansions(), SameTypeOfBoundaryConditions(), SetUpDG(), v_EvaluateBoundaryConditions(), Nektar::MultiRegions::ContField3D::v_FillBndCondFromField(), v_GetBndCondExpansions(), v_GetBoundaryToElmtMap(), v_GetFwdBwdTracePhys(), v_GetRobinBCInfo(), v_HelmSolve(), Nektar::MultiRegions::ContField3D::v_HelmSolve(), Nektar::MultiRegions::ContField3D::v_ImposeDirichletConditions(), and v_UpdateBndCondExpansion().

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

An array which contains the information about the boundary condition on the different boundary regions.

Definition at line 98 of file DisContField3D.h.

Referenced by Nektar::MultiRegions::ContField3D::ContField3D(), GenerateBoundaryConditionExpansion(), Nektar::MultiRegions::ContField3D::GenerateDirBndCondForcing(), SameTypeOfBoundaryConditions(), SetUpDG(), v_EvaluateBoundaryConditions(), v_GetBndConditions(), v_GetBoundaryToElmtMap(), v_GetFwdBwdTracePhys(), v_GetRobinBCInfo(), v_HelmSolve(), Nektar::MultiRegions::ContField3D::v_HelmSolve(), Nektar::MultiRegions::ContField3D::v_ImposeDirichletConditions(), and v_UpdateBndConditions().

std::set<int> Nektar::MultiRegions::DisContField3D::m_boundaryFaces
protected

A set storing the global IDs of any boundary faces.

Definition at line 107 of file DisContField3D.h.

Referenced by IsLeftAdjacentFace(), and SetUpDG().

GlobalLinSysMapShPtr Nektar::MultiRegions::DisContField3D::m_globalBndMat
protected

Definition at line 100 of file DisContField3D.h.

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

vector<bool> Nektar::MultiRegions::DisContField3D::m_leftAdjacentFaces
protected

Definition at line 128 of file DisContField3D.h.

Referenced by SetUpDG(), and v_GetFwdBwdTracePhys().

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

Definition at line 136 of file DisContField3D.h.

Referenced by SetUpDG(), and v_GetFwdBwdTracePhys().

PeriodicMap Nektar::MultiRegions::DisContField3D::m_periodicEdges
protected

A map which identifies groups of periodic edges.

Definition at line 117 of file DisContField3D.h.

Referenced by Nektar::MultiRegions::ContField3D::ContField3D(), DisContField3D(), FindPeriodicFaces(), and v_GetPeriodicEntities().

PeriodicMap Nektar::MultiRegions::DisContField3D::m_periodicFaces
protected

A map which identifies pairs of periodic faces.

Definition at line 112 of file DisContField3D.h.

Referenced by Nektar::MultiRegions::ContField3D::ContField3D(), DisContField3D(), FindPeriodicFaces(), IsLeftAdjacentFace(), SetUpDG(), and v_GetPeriodicEntities().

vector<int> Nektar::MultiRegions::DisContField3D::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 135 of file DisContField3D.h.

Referenced by SetUpDG(), and v_GetFwdBwdTracePhys().

PeriodicMap Nektar::MultiRegions::DisContField3D::m_periodicVerts
protected

A map which identifies groups of periodic vertices.

Definition at line 122 of file DisContField3D.h.

Referenced by Nektar::MultiRegions::ContField3D::ContField3D(), DisContField3D(), FindPeriodicFaces(), and v_GetPeriodicEntities().

ExpListSharedPtr Nektar::MultiRegions::DisContField3D::m_trace
protected

Definition at line 101 of file DisContField3D.h.

Referenced by DisContField3D(), EvaluateHDGPostProcessing(), IsLeftAdjacentFace(), SetUpDG(), v_AddFwdBwdTraceIntegral(), v_AddTraceIntegral(), v_ExtractTracePhys(), v_GetFwdBwdTracePhys(), v_GetTrace(), and v_HelmSolve().

AssemblyMapDGSharedPtr Nektar::MultiRegions::DisContField3D::m_traceMap
protected

Definition at line 102 of file DisContField3D.h.

Referenced by DisContField3D(), EvaluateHDGPostProcessing(), GetGlobalBndLinSys(), IsLeftAdjacentFace(), SetUpDG(), v_AddFwdBwdTraceIntegral(), v_AddTraceIntegral(), v_ExtractTracePhys(), v_GeneralMatrixOp(), v_GetFwdBwdTracePhys(), v_GetTraceMap(), and v_HelmSolve().