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

This class is the abstraction of a global discontinuous two- dimensional spectral/hp element expansion which approximates the solution of a set of partial differential equations. More...

#include <DisContField1D.h>

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

Public Member Functions

 DisContField1D ()
 Default constructor.
 DisContField1D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &graph1D, const std::string &variable, const bool SetUpJustDG=true)
 Constructs a 1D discontinuous field based on a mesh and boundary conditions.
 DisContField1D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &graph1D, const SpatialDomains::CompositeMap &domain, const SpatialDomains::BoundaryConditions &Allbcs, const std::string &variable, bool SetToOneSpaceDimensions=false)
 Constructor for a DisContField1D from a List of subdomains New Constructor for arterial network.
 DisContField1D (const DisContField1D &In)
 Constructs a 1D discontinuous field based on an existing field.
 DisContField1D (const ExpList1D &In)
 Constructs a 1D discontinuous field based on an existing field. (needed in order to use ContField( const ExpList1D &In) constructor.
virtual ~DisContField1D ()
 Destructor.
GlobalLinSysSharedPtr GetGlobalBndLinSys (const GlobalLinSysKey &mkey)
 For a given key, returns the associated global linear system.
vector< bool > & GetNegatedFluxNormal (void)
- Public Member Functions inherited from Nektar::MultiRegions::ExpList1D
 ExpList1D ()
 The default constructor.
 ExpList1D (const ExpList1D &In, const bool DeclareCoeffPhysArrays=true)
 The copy constructor.
 ExpList1D (const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::BasisKey &Ba, const SpatialDomains::MeshGraphSharedPtr &graph1D)
 Construct an ExpList1D from a given graph.
 ExpList1D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &graph1D, const bool DeclareCoeffPhysArrays=true)
 This constructor sets up a list of local expansions based on an input graph1D.
 ExpList1D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &graph1D, const SpatialDomains::CompositeMap &domain, const bool DeclareCoeffPhysArrays=true, const std::string var="DefaultVar", bool SetToOneSpaceDimension=false)
 This constructor sets up a list of local expansions based on an input compositeMap.
 ExpList1D (const SpatialDomains::CompositeMap &domain, const SpatialDomains::MeshGraphSharedPtr &graph2D, const bool DeclareCoeffPhysArrays=true, const std::string variable="DefaultVar")
 Specialised constructor for Neumann boundary conditions in DisContField2D and ContField2D.
 ExpList1D (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::CompositeMap &domain, const SpatialDomains::MeshGraphSharedPtr &graph1D, int i, const bool DeclareCoeffPhysArrays=true)
 ExpList1D (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, const ExpListSharedPtr > &bndConstraint, const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > &bndCond, const LocalRegions::ExpansionVector &locexp, const SpatialDomains::MeshGraphSharedPtr &graph2D, const PeriodicMap &periodicEdges, const bool DeclareCoeffPhysArrays=true, const std::string variable="DefaultVar")
 Specialised constructor for trace expansions.
virtual ~ExpList1D ()
 Destructor.
void PostProcess (LibUtilities::KernelSharedPtr kernel, Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray, NekDouble h, int elmId=0)
 Performs the post-processing on a specified element.
void PeriodicEval (Array< OneD, NekDouble > &inarray1, Array< OneD, NekDouble > &inarray2, NekDouble h, int nmodes, Array< OneD, NekDouble > &outarray)
 Evaluates the global spectral/hp expansion at some arbitray set of points.
void ParNormalSign (Array< OneD, NekDouble > &normsign)
 Set up the normals on each expansion.
- Public Member Functions inherited from Nektar::MultiRegions::ExpList
 ExpList ()
 The default constructor.
 ExpList (const LibUtilities::SessionReaderSharedPtr &pSession)
 The default constructor.
 ExpList (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 The default constructor.
 ExpList (const ExpList &in, const bool DeclareCoeffPhysArrays=true)
 The copy constructor.
virtual ~ExpList ()
 The default destructor.
int GetNcoeffs (void) const
 Returns the total number of local degrees of freedom $N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_m$.
int GetNcoeffs (const int eid) const
 Returns the total number of local degrees of freedom for element eid.
ExpansionType GetExpType (void)
 Returns the type of the expansion.
void SetExpType (ExpansionType Type)
 Returns the type of the expansion.
int EvalBasisNumModesMax (void) const
 Evaulates the maximum number of modes in the elemental basis order over all elements.
const Array< OneD, int > EvalBasisNumModesMaxPerExp (void) const
 Returns the vector of the number of modes in the elemental basis order over all elements.
int GetTotPoints (void) const
 Returns the total number of quadrature points m_npoints $=Q_{\mathrm{tot}}$.
int GetTotPoints (const int eid) const
 Returns the total number of quadrature points for eid's element $=Q_{\mathrm{tot}}$.
int GetNpoints (void) const
 Returns the total number of quadrature points m_npoints $=Q_{\mathrm{tot}}$.
int Get1DScaledTotPoints (const NekDouble scale) const
 Returns the total number of qudature points scaled by the factor scale on each 1D direction.
void SetWaveSpace (const bool wavespace)
 Sets the wave space to the one of the possible configuration true or false.
void SetModifiedBasis (const bool modbasis)
 Set Modified Basis for the stability analysis.
void SetPhys (int i, NekDouble val)
 Set the i th value of m_phys to value val.
bool GetWaveSpace (void) const
 This function returns the third direction expansion condition, which can be in wave space (coefficient) or not It is stored in the variable m_WaveSpace.
void SetPhys (const Array< OneD, const NekDouble > &inarray)
 Fills the array m_phys.
void SetPhysArray (Array< OneD, NekDouble > &inarray)
 Sets the array m_phys.
void SetPhysState (const bool physState)
 This function manually sets whether the array of physical values $\boldsymbol{u}_l$ (implemented as m_phys) is filled or not.
bool GetPhysState (void) const
 This function indicates whether the array of physical values $\boldsymbol{u}_l$ (implemented as m_phys) is filled or not.
NekDouble PhysIntegral (void)
 This function integrates a function $f(\boldsymbol{x})$ over the domain consisting of all the elements of the expansion.
NekDouble PhysIntegral (const Array< OneD, const NekDouble > &inarray)
 This function integrates a function $f(\boldsymbol{x})$ over the domain consisting of all the elements of the expansion.
void IProductWRTBase_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function calculates the inner product of a function $f(\boldsymbol{x})$ with respect to all {local} expansion modes $\phi_n^e(\boldsymbol{x})$.
void IProductWRTBase (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
void IProductWRTDerivBase (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function calculates the inner product of a function $f(\boldsymbol{x})$ with respect to the derivative (in direction.
void FwdTrans_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function elementally evaluates the forward transformation of a function $u(\boldsymbol{x})$ onto the global spectral/hp expansion.
void FwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
void MultiplyByElmtInvMass (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function elementally mulplies the coefficient space of Sin my the elemental inverse of the mass matrix.
void MultiplyByInvMassMatrix (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
void SmoothField (Array< OneD, NekDouble > &field)
 Smooth a field across elements.
void HelmSolve (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const FlagList &flags, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff=StdRegions::NullVarCoeffMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 Solve helmholtz problem.
void LinearAdvectionDiffusionReactionSolve (const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 Solve Advection Diffusion Reaction.
void LinearAdvectionReactionSolve (const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
 Solve Advection Diffusion Reaction.
void FwdTrans_BndConstrained (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void BwdTrans_IterPerExp (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function elementally evaluates the backward transformation of the global spectral/hp element expansion.
void BwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
void GetCoords (Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1=NullNekDouble1DArray, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
 This function calculates the coordinates of all the elemental quadrature points $\boldsymbol{x}_i$.
void HomogeneousFwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
void HomogeneousBwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
void DealiasedProd (const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
void GetBCValues (Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
void NormVectorIProductWRTBase (Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
void ApplyGeomInfo ()
 Apply geometry information to each expansion.
void WriteTecplotHeader (std::ofstream &outfile, std::string var="")
void WriteTecplotZone (std::ofstream &outfile, int expansion=-1)
void WriteTecplotField (std::ofstream &outfile, int expansion=-1)
void WriteTecplotConnectivity (std::ofstream &outfile, int expansion=-1)
void WriteVtkHeader (std::ofstream &outfile)
void WriteVtkFooter (std::ofstream &outfile)
void WriteVtkPieceHeader (std::ofstream &outfile, int expansion)
void WriteVtkPieceFooter (std::ofstream &outfile, int expansion)
void WriteVtkPieceData (std::ofstream &outfile, int expansion, std::string var="v")
int GetCoordim (int eid)
 This function returns the dimension of the coordinates of the element eid.
void SetCoeff (int i, NekDouble val)
 Set the i th coefficiient in m_coeffs to value val.
void SetCoeffs (int i, NekDouble val)
 Set the i th coefficiient in m_coeffs to value val.
void SetCoeffsArray (Array< OneD, NekDouble > &inarray)
 Set the m_coeffs array to inarray.
const Array< OneD, const
NekDouble > & 
GetCoeffs () const
 This function returns (a reference to) the array $\boldsymbol{\hat{u}}_l$ (implemented as m_coeffs) containing all local expansion coefficients.
void ImposeDirichletConditions (Array< OneD, NekDouble > &outarray)
 Impose Dirichlet Boundary Conditions onto Array.
void FillBndCondFromField (void)
 Fill Bnd Condition expansion from the values stored in expansion.
void LocalToGlobal (void)
 Put the coefficients into global ordering using m_coeffs.
void GlobalToLocal (void)
 Put the coefficients into local ordering and place in m_coeffs.
NekDouble GetCoeff (int i)
 Get the i th value (coefficient) of m_coeffs.
NekDouble GetCoeffs (int i)
 Get the i th value (coefficient) of m_coeffs.
const Array< OneD, const
NekDouble > & 
GetPhys () const
 This function returns (a reference to) the array $\boldsymbol{u}_l$ (implemented as m_phys) containing the function $u^{\delta}(\boldsymbol{x})$ evaluated at the quadrature points.
NekDouble Linf (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
 This function calculates the $L_\infty$ error of the global spectral/hp element approximation.
NekDouble L2 (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
 This function calculates the $L_2$ error with respect to soln of the global spectral/hp element approximation.
NekDouble H1 (const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
 Calculates the $H^1$ error of the global spectral/hp element approximation.
NekDouble Integral (const Array< OneD, const NekDouble > &inarray)
Array< OneD, const NekDoubleHomogeneousEnergy (void)
 This function calculates the energy associated with each one of the modesof a 3D homogeneous nD expansion.
void SetHomo1DSpecVanVisc (Array< OneD, NekDouble > visc)
 This function sets the Spectral Vanishing Viscosity in homogeneous1D expansion.
Array< OneD, const unsigned int > GetZIDs (void)
 This function returns a vector containing the wave numbers in z-direction associated with the 3D homogenous expansion. Required if a parellelisation is applied in the Fourier direction.
LibUtilities::TranspositionSharedPtr GetTransposition (void)
 This function returns the transposition class associaed with the homogeneous expansion.
NekDouble GetHomoLen (void)
 This function returns the Width of homogeneous direction associaed with the homogeneous expansion.
Array< OneD, const unsigned int > GetYIDs (void)
 This function returns a vector containing the wave numbers in y-direction associated with the 3D homogenous expansion. Required if a parellelisation is applied in the Fourier direction.
void PhysInterp1DScaled (const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function interpolates the physical space points in inarray to outarray using the same points defined in the expansion but where the number of points are rescaled by 1DScale.
void PhysGalerkinProjection1DScaled (const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This function Galerkin projects the physical space points in inarray to outarray where inarray is assumed to be defined in the expansion but where the number of points are rescaled by 1DScale.
int GetExpSize (void)
 This function returns the number of elements in the expansion.
int GetNumElmts (void)
 This function returns the number of elements in the expansion which may be different for a homogeoenous extended expansionp.
const boost::shared_ptr
< LocalRegions::ExpansionVector
GetExp () const
 This function returns the vector of elements in the expansion.
LocalRegions::ExpansionSharedPtrGetExp (int n) const
 This function returns (a shared pointer to) the local elemental expansion of the $n^{\mathrm{th}}$ element.
LocalRegions::ExpansionSharedPtrGetExp (const Array< OneD, const NekDouble > &gloCoord)
 This function returns (a shared pointer to) the local elemental expansion containing the arbitrary point given by gloCoord.
int GetExpIndex (const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0, bool returnNearestElmt=false)
int GetExpIndex (const Array< OneD, const NekDouble > &gloCoords, Array< OneD, NekDouble > &locCoords, NekDouble tol=0.0, bool returnNearestElmt=false)
int GetCoeff_Offset (int n) const
 Get the start offset position for a global list of m_coeffs correspoinding to element n.
int GetPhys_Offset (int n) const
 Get the start offset position for a global list of m_phys correspoinding to element n.
int GetOffset_Elmt_Id (int n) const
 Get the element id associated with the n th consecutive block of data in m_phys and m_coeffs.
Array< OneD, NekDouble > & UpdateCoeffs ()
 This function returns (a reference to) the array $\boldsymbol{\hat{u}}_l$ (implemented as m_coeffs) containing all local expansion coefficients.
Array< OneD, NekDouble > & UpdatePhys ()
 This function returns (a reference to) the array $\boldsymbol{u}_l$ (implemented as m_phys) containing the function $u^{\delta}(\boldsymbol{x})$ evaluated at the quadrature points.
void PhysDeriv (Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
void PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1=NullNekDouble1DArray, Array< OneD, NekDouble > &out_d2=NullNekDouble1DArray)
 This function discretely evaluates the derivative of a function $f(\boldsymbol{x})$ on the domain consisting of all elements of the expansion.
void PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
const Array< OneD, const
boost::shared_ptr< ExpList > > & 
GetBndCondExpansions ()
boost::shared_ptr< ExpList > & UpdateBndCondExpansion (int i)
void Upwind (const Array< OneD, const Array< OneD, NekDouble > > &Vec, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
void Upwind (const Array< OneD, const NekDouble > &Vn, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
boost::shared_ptr< ExpList > & GetTrace ()
boost::shared_ptr
< AssemblyMapDG > & 
GetTraceMap (void)
const Array< OneD, const int > & GetTraceBndMap (void)
void GetNormals (Array< OneD, Array< OneD, NekDouble > > &normals)
void AddTraceIntegral (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
void AddTraceIntegral (const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
void AddFwdBwdTraceIntegral (const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
void GetFwdBwdTracePhys (Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
void GetFwdBwdTracePhys (const Array< OneD, const NekDouble > &field, Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
void ExtractTracePhys (Array< OneD, NekDouble > &outarray)
void ExtractTracePhys (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
const Array< OneD, const
SpatialDomains::BoundaryConditionShPtr > & 
GetBndConditions ()
Array< OneD,
SpatialDomains::BoundaryConditionShPtr > & 
UpdateBndConditions ()
void EvaluateBoundaryConditions (const NekDouble time=0.0, const std::string varName="", const NekDouble=NekConstants::kNekUnsetDouble, const NekDouble=NekConstants::kNekUnsetDouble)
void GeneralMatrixOp (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
 This function calculates the result of the multiplication of a matrix of type specified by mkey with a vector given by inarray.
void GeneralMatrixOp_IterPerExp (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void SetUpPhysNormals ()
void GetBoundaryToElmtMap (Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
void GeneralGetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef, int NumHomoDir=0, Array< OneD, LibUtilities::BasisSharedPtr > &HomoBasis=LibUtilities::NullBasisSharedPtr1DArray, std::vector< NekDouble > &HomoLen=LibUtilities::NullNekDoubleVector, std::vector< unsigned int > &HomoZIDs=LibUtilities::NullUnsignedIntVector, std::vector< unsigned int > &HomoYIDs=LibUtilities::NullUnsignedIntVector)
const
NekOptimize::GlobalOptParamSharedPtr
GetGlobalOptParam (void)
map< int, RobinBCInfoSharedPtrGetRobinBCInfo ()
void GetPeriodicEntities (PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces=NullPeriodicMap)
std::vector
< LibUtilities::FieldDefinitionsSharedPtr
GetFieldDefinitions ()
void GetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
void AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
 Append the element data listed in elements fielddef->m_ElementIDs onto fielddata.
void AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, Array< OneD, NekDouble > &coeffs)
 Append the data in coeffs listed in elements fielddef->m_ElementIDs onto fielddata.
void ExtractElmtDataToCoeffs (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
 Extract the data in fielddata into the coeffs using the basic ExpList Elemental expansions rather than planes in homogeneous case.
void ExtractCoeffsToCoeffs (const boost::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
 Extract the data from fromField using fromExpList the coeffs using the basic ExpList Elemental expansions rather than planes in homogeneous case.
void ExtractDataToCoeffs (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
 Extract the data in fielddata into the coeffs.
boost::shared_ptr< ExpListGetSharedThisPtr ()
 Returns a shared pointer to the current object.
boost::shared_ptr
< LibUtilities::SessionReader
GetSession ()
 Returns the session object.
boost::shared_ptr
< LibUtilities::Comm
GetComm ()
 Returns the comm object.
SpatialDomains::MeshGraphSharedPtr GetGraph ()
LibUtilities::BasisSharedPtr GetHomogeneousBasis (void)
boost::shared_ptr< ExpList > & GetPlane (int n)

Protected Member Functions

void GenerateBoundaryConditionExpansion (const SpatialDomains::MeshGraphSharedPtr &graph1D, const SpatialDomains::BoundaryConditions &bcs, const std::string variable)
 Discretises the boundary conditions.
void FindPeriodicVertices (const SpatialDomains::BoundaryConditions &bcs, const std::string variable)
 Generate a associative map of periodic vertices in a mesh.
virtual ExpListSharedPtrv_GetTrace ()
virtual AssemblyMapDGSharedPtrv_GetTraceMap (void)
virtual void v_AddTraceIntegral (const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
virtual void v_GetFwdBwdTracePhys (Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
virtual void v_GetFwdBwdTracePhys (const Array< OneD, const NekDouble > &field, Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
 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_ExtractTracePhys (Array< OneD, NekDouble > &outarray)
virtual void v_ExtractTracePhys (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
 This method extracts the trace (verts in 1D) from the field inarray and puts the values in outarray.
void SetBoundaryConditionExpansion (const SpatialDomains::MeshGraphSharedPtr &graph1D, const SpatialDomains::BoundaryConditions &bcs, const std::string variable, Array< OneD, MultiRegions::ExpListSharedPtr > &bndCondExpansions, Array< OneD, SpatialDomains::BoundaryConditionShPtr > &bndConditions)
 Populates the list of boundary condition expansions.
void SetMultiDomainBoundaryConditionExpansion (const SpatialDomains::MeshGraphSharedPtr &graph1D, const SpatialDomains::BoundaryConditions &bcs, const std::string variable, Array< OneD, MultiRegions::ExpListSharedPtr > &bndCondExpansions, Array< OneD, SpatialDomains::BoundaryConditionShPtr > &bndConditions, int subdomain)
 Populates the list of boundary condition expansions in multidomain case.
void GenerateFieldBnd1D (SpatialDomains::BoundaryConditions &bcs, const std::string variable)
virtual map< int,
RobinBCInfoSharedPtr
v_GetRobinBCInfo ()
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_GetBoundaryToElmtMap (Array< OneD, int > &ElmtID, Array< OneD, int > &VertID)
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)
 Evaluate all boundary conditions at a given time..
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)
 Solve the Helmholtz equation.
- Protected Member Functions inherited from Nektar::MultiRegions::ExpList1D
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)
 Upwind the Fwd and Bwd states based on the velocity field given by Vec.
void v_Upwind (const Array< OneD, const NekDouble > &Vn, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
 Upwind the Fwd and Bwd states based on the one- dimensional normal velocity field given by Vn.
void v_GetNormals (Array< OneD, Array< OneD, NekDouble > > &normals)
 Populate normals with the normals of all expansions.
- 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 const Array< OneD,
const int > & 
v_GetTraceBndMap ()
virtual void v_AddTraceIntegral (const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
virtual void v_AddFwdBwdTraceIntegral (const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
virtual void v_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_GeneralMatrixOp (const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
virtual void v_GetCoords (Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
virtual void v_PhysDeriv (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
virtual void v_PhysDeriv (const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
virtual void v_PhysDeriv (Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
virtual void v_HomogeneousFwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
virtual void v_HomogeneousBwdTrans (const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal, bool Shuff=true, bool UnShuff=true)
virtual void v_DealiasedProd (const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
virtual void v_GetBCValues (Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
virtual void v_NormVectorIProductWRTBase (Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
virtual void v_SetUpPhysNormals ()
virtual void v_ReadGlobalOptimizationParameters ()
virtual std::vector
< LibUtilities::FieldDefinitionsSharedPtr
v_GetFieldDefinitions (void)
virtual void v_GetFieldDefinitions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
virtual void v_AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
virtual void v_AppendFieldData (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, Array< OneD, NekDouble > &coeffs)
virtual void v_ExtractDataToCoeffs (LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
 Extract data from raw field data into expansion list.
virtual void v_ExtractCoeffsToCoeffs (const boost::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
virtual void v_WriteTecplotHeader (std::ofstream &outfile, std::string var="")
virtual void v_WriteTecplotZone (std::ofstream &outfile, int expansion)
virtual void v_WriteTecplotField (std::ofstream &outfile, int expansion)
virtual void v_WriteTecplotConnectivity (std::ofstream &outfile, int expansion)
virtual void v_WriteVtkPieceHeader (std::ofstream &outfile, int expansion)
virtual void v_WriteVtkPieceData (std::ofstream &outfile, int expansion, std::string var)
virtual NekDouble v_L2 (const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
virtual NekDouble v_Integral (const Array< OneD, const NekDouble > &inarray)
virtual Array< OneD, const
NekDouble
v_HomogeneousEnergy (void)
virtual
LibUtilities::TranspositionSharedPtr 
v_GetTransposition (void)
virtual NekDouble v_GetHomoLen (void)
virtual Array< OneD, const
unsigned int > 
v_GetZIDs (void)
virtual Array< OneD, const
unsigned int > 
v_GetYIDs (void)
virtual void v_PhysInterp1DScaled (const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_PhysGalerkinProjection1DScaled (const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void ExtractFileBCs (const std::string &fileName, const std::string &varName, const boost::shared_ptr< ExpList > locExpList)

Protected Attributes

int m_numDirBndCondExpansions
 The number of boundary segments on which Dirichlet boundary conditions are imposed.
Array< OneD,
MultiRegions::ExpListSharedPtr
m_bndCondExpansions
 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
 Global boundary matrix.
ExpListSharedPtr m_trace
 Trace space storage for points between elements.
AssemblyMapDGSharedPtr m_traceMap
 Local to global DG mapping for trace space.
std::set< int > m_boundaryVerts
 A set storing the global IDs of any boundary edges.
PeriodicMap m_periodicVerts
 A map which identifies groups of periodic vertices.
vector< int > m_periodicFwdCopy
 A vector indicating degress of freedom which need to be copied from forwards to backwards space in case of a periodic boundary condition.
vector< int > m_periodicBwdCopy
vector< bool > m_leftAdjacentVerts

Private Member Functions

void SetUpDG (const std::string &variable)
bool IsLeftAdjacentVertex (const int n, const int e)
SpatialDomains::BoundaryConditionsSharedPtr GetDomainBCs (const SpatialDomains::CompositeMap &domain, const SpatialDomains::BoundaryConditions &Allbcs, const std::string &variable)

Private Attributes

vector< bool > m_negatedFluxNormal

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

This class is the abstraction of a global discontinuous two- dimensional spectral/hp element expansion which approximates the solution of a set of partial differential equations.

This class augments the list of local expansions inherited from ExpList1D with boundary conditions. Inter-element boundaries are handled using an discontinuous Galerkin scheme.

Definition at line 55 of file DisContField1D.h.

Constructor & Destructor Documentation

Nektar::MultiRegions::DisContField1D::DisContField1D ( )

Default constructor.

Constructs an empty expansion list with no boundary conditions.

Definition at line 56 of file DisContField1D.cpp.

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

Constructs a 1D discontinuous field based on a mesh and boundary conditions.

An expansion list for the boundary expansions is generated first for the field. These are subsequently evaluated for time zero. The trace map is then constructed.

Parameters
graph1DA mesh, containing information about the domain and the spectral/hp element expansions.
bcsInformation about the enforced boundary conditions.
variableThe session variable associated with the boundary conditions to enforce.
solnTypeType of global system to use.

Definition at line 76 of file DisContField1D.cpp.

References Nektar::MultiRegions::ExpList::ApplyGeomInfo(), Nektar::MultiRegions::ExpList::EvaluateBoundaryConditions(), FindPeriodicVertices(), GenerateBoundaryConditionExpansion(), Nektar::MultiRegions::ExpList::m_session, and SetUpDG().

: ExpList1D(pSession,graph1D),
{
GenerateBoundaryConditionExpansion(graph1D,bcs,variable);
FindPeriodicVertices(bcs,variable);
if(SetUpJustDG)
{
SetUpDG(variable);
}
}
Nektar::MultiRegions::DisContField1D::DisContField1D ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr graph1D,
const SpatialDomains::CompositeMap domain,
const SpatialDomains::BoundaryConditions Allbcs,
const std::string &  variable,
bool  SetToOneSpaceDimension = false 
)

Constructor for a DisContField1D from a List of subdomains New Constructor for arterial network.

Constructor for use in multidomain computations where a domain list can be passed instead of graph1D

Parameters
domainSubdomain specified in the inputfile from which the DisContField1D is set up

Definition at line 418 of file DisContField1D.cpp.

References Nektar::MultiRegions::ExpList::ApplyGeomInfo(), Nektar::MultiRegions::ExpList::EvaluateBoundaryConditions(), FindPeriodicVertices(), GenerateBoundaryConditionExpansion(), GetDomainBCs(), and SetUpDG().

:
ExpList1D(pSession,graph1D,domain, true,variable,SetToOneSpaceDimension),
{
SpatialDomains::BoundaryConditionsSharedPtr DomBCs = GetDomainBCs(domain,Allbcs,variable);
GenerateBoundaryConditionExpansion(graph1D,*DomBCs,variable);
FindPeriodicVertices(*DomBCs,variable);
SetUpDG(variable);
}
Nektar::MultiRegions::DisContField1D::DisContField1D ( const DisContField1D In)

Constructs a 1D discontinuous field based on an existing field.

Constructs a field as a copy of an existing field.

Parameters
InExisting DisContField1D object to copy.

Definition at line 443 of file DisContField1D.cpp.

:
ExpList1D(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_boundaryVerts(In.m_boundaryVerts),
m_periodicVerts(In.m_periodicVerts),
m_periodicFwdCopy(In.m_periodicFwdCopy),
m_periodicBwdCopy(In.m_periodicBwdCopy),
m_leftAdjacentVerts(In.m_leftAdjacentVerts)
{
}
Nektar::MultiRegions::DisContField1D::DisContField1D ( const ExpList1D In)

Constructs a 1D discontinuous field based on an existing field. (needed in order to use ContField( const ExpList1D &In) constructor.

Constructs a field as a copy of an existing explist1D field.

Parameters
InExisting ExpList1D object to copy.

Definition at line 463 of file DisContField1D.cpp.

:
{
}
Nektar::MultiRegions::DisContField1D::~DisContField1D ( )
virtual

Destructor.

Definition at line 471 of file DisContField1D.cpp.

{
}

Member Function Documentation

void Nektar::MultiRegions::DisContField1D::FindPeriodicVertices ( const SpatialDomains::BoundaryConditions bcs,
const std::string  variable 
)
protected

Generate a associative map of periodic vertices in a mesh.

Parameters
graph1DA mesh containing information about the domain and the spectral/hp element expansion.
bcsInformation about the boundary conditions.
variableSpecifies the field.
periodicVerticesMap into which the list of periodic vertices is placed.

Definition at line 535 of file DisContField1D.cpp.

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

Referenced by DisContField1D().

{
= boost::dynamic_pointer_cast<
SpatialDomains::BoundaryRegionCollection::const_iterator it;
m_session->GetComm()->GetRowComm();
int i, region1ID, region2ID;
map<int,int> BregionToVertMap;
// Construct list of all periodic Region and their global vertex on
// this process.
for (it = bregions.begin(); it != bregions.end(); ++it)
{
locBCond = GetBoundaryCondition(bconditions, it->first, variable);
if (locBCond->GetBoundaryConditionType()
{
continue;
}
int id = (*(it->second->begin()->second))[0]->GetGlobalID();
BregionToVertMap[it->first] = id;
}
set<int> islocal;
int n = vComm->GetSize();
int p = vComm->GetRank();
Array<OneD, int> nregions(n, 0);
nregions[p] = BregionToVertMap.size();
vComm->AllReduce(nregions, LibUtilities::ReduceSum);
int totRegions = Vmath::Vsum(n, nregions, 1);
Array<OneD, int> regOffset(n, 0);
for (i = 1; i < n; ++i)
{
regOffset[i] = regOffset[i-1] + nregions[i-1];
}
Array<OneD, int> bregmap(totRegions, 0);
Array<OneD, int> bregid (totRegions, 0);
for(i = regOffset[p], iit = BregionToVertMap.begin();
iit != BregionToVertMap.end(); ++iit, ++i)
{
bregid [i] = iit->first;
bregmap[i] = iit->second;
islocal.insert(iit->first);
}
vComm->AllReduce(bregmap, LibUtilities::ReduceSum);
vComm->AllReduce(bregid, LibUtilities::ReduceSum);
for (int i = 0; i < totRegions; ++i)
{
BregionToVertMap[bregid[i]] = bregmap[i];
}
// Construct list of all periodic pairs local to this process.
for (it = bregions.begin(); it != bregions.end(); ++it)
{
locBCond = GetBoundaryCondition(bconditions, it->first, variable);
if (locBCond->GetBoundaryConditionType()
{
continue;
}
// Identify periodic boundary region IDs.
region1ID = it->first;
region2ID = boost::static_pointer_cast<
locBCond)->m_connectedBoundaryRegion;
ASSERTL0(BregionToVertMap.count(region1ID) != 0,
"Cannot determine vertex of region1ID");
ASSERTL0(BregionToVertMap.count(region2ID) != 0,
"Cannot determine vertex of region2ID");
PeriodicEntity ent(BregionToVertMap[region2ID],
islocal.count(region2ID) != 0);
m_periodicVerts[BregionToVertMap[region1ID]].push_back(ent);
}
}
void Nektar::MultiRegions::DisContField1D::GenerateBoundaryConditionExpansion ( const SpatialDomains::MeshGraphSharedPtr graph1D,
const SpatialDomains::BoundaryConditions bcs,
const std::string  variable 
)
protected

Discretises the boundary conditions.

Generate the boundary condition expansion list

Parameters
graph1DA mesh, containing information about the domain and the spectral/hp element expansions.
bcsInformation about the enforced boundary conditions.
variableThe session variable associated with the boundary conditions to enforce.

Definition at line 485 of file DisContField1D.cpp.

References Nektar::SpatialDomains::ePeriodic, Nektar::MultiRegions::ExpList::GetBoundaryCondition(), Nektar::SpatialDomains::BoundaryConditions::GetBoundaryConditions(), Nektar::SpatialDomains::BoundaryConditions::GetBoundaryRegions(), Nektar::iterator, m_bndCondExpansions, m_bndConditions, and SetBoundaryConditionExpansion().

Referenced by DisContField1D().

{
int cnt = 0;
SpatialDomains::BoundaryRegionCollection::const_iterator it;
// count the number of non-periodic boundary points
for (it = bregions.begin(); it != bregions.end(); ++it)
{
const SpatialDomains::BoundaryConditionShPtr boundaryCondition =
GetBoundaryCondition(bconditions, it->first, variable);
if (boundaryCondition->GetBoundaryConditionType() !=
{
for (bregionIt = it->second->begin();
bregionIt != it->second->end(); bregionIt++)
{
cnt += bregionIt->second->size();
}
}
}
= Array<OneD,MultiRegions::ExpListSharedPtr>(cnt);
= Array<OneD,SpatialDomains::BoundaryConditionShPtr>(cnt);
SetBoundaryConditionExpansion(graph1D,bcs,variable,
m_bndConditions);
}
void Nektar::MultiRegions::DisContField1D::GenerateFieldBnd1D ( SpatialDomains::BoundaryConditions bcs,
const std::string  variable 
)
protected
SpatialDomains::BoundaryConditionsSharedPtr Nektar::MultiRegions::DisContField1D::GetDomainBCs ( const SpatialDomains::CompositeMap domain,
const SpatialDomains::BoundaryConditions Allbcs,
const std::string &  variable 
)
private

Definition at line 308 of file DisContField1D.cpp.

References ASSERTL1, Nektar::SpatialDomains::BoundaryConditions::GetBoundaryConditions(), Nektar::SpatialDomains::BoundaryConditions::GetBoundaryRegions(), Nektar::iterator, and Nektar::MultiRegions::ExpList::m_session.

Referenced by DisContField1D().

{
map<int,int> GeometryToRegionsMap;
SpatialDomains::BoundaryRegionCollection::const_iterator it;
= Allbcs.GetBoundaryRegions();
// Set up a map of all boundary regions
for(it = bregions.begin(); it != bregions.end(); ++it)
{
for (bregionIt = it->second->begin();
bregionIt != it->second->end(); bregionIt++)
{
// can assume that all regions only contain one point in 1D
// Really do not need loop above
int id = (*(bregionIt->second))[0]->GetGlobalID();
GeometryToRegionsMap[id] = it->first;
}
}
map<int,SpatialDomains::GeometrySharedPtr> EndOfDomain;
// Now find out which points in domain have only one vertex
for(domIt = domain.begin(); domIt != domain.end(); ++domIt)
{
SpatialDomains::Composite geomvector = domIt->second;
for(int i = 0; i < geomvector->size(); ++i)
{
for(int j = 0; j < 2; ++j)
{
int vid = (*geomvector)[i]->GetVid(j);
if(EndOfDomain.count(vid) == 0)
{
EndOfDomain[vid] = (*geomvector)[i]->GetVertex(j);
}
else
{
EndOfDomain.erase(vid);
}
}
}
}
ASSERTL1(EndOfDomain.size() == 2,"Did not find two ends of domain");
int numNewBc = 1;
for(regIt = EndOfDomain.begin(); regIt != EndOfDomain.end(); ++regIt)
{
if(GeometryToRegionsMap.count(regIt->first) != 0) // Set up boundary condition up
{
map<int,int>::iterator iter = GeometryToRegionsMap.find(regIt->first);
ASSERTL1(iter != GeometryToRegionsMap.end(),"Failied to find GeometryToRegionMap");
int regionId = iter->second;
SpatialDomains::BoundaryRegionCollection::const_iterator bregionsIter = bregions.find(regionId);
ASSERTL1(bregionsIter != bregions.end(),"Failed to find boundary region");
SpatialDomains::BoundaryRegionShPtr breg = bregionsIter->second;
returnval->AddBoundaryRegions (regionId,breg);
SpatialDomains::BoundaryConditionCollection::const_iterator bconditionsIter = bconditions.find(regionId);
ASSERTL1(bconditionsIter != bconditions.end(),"Failed to find boundary collection");
SpatialDomains::BoundaryConditionMapShPtr bcond = bconditionsIter->second;
returnval->AddBoundaryConditions(regionId,bcond);
}
else // Set up an undefined region.
{
// Set up Composite (GemetryVector) to contain vertex and put into bRegion
gvec->push_back(regIt->second);
(*breg)[regIt->first] = gvec;
returnval->AddBoundaryRegions(bregions.size()+numNewBc,breg);
// Set up just boundary condition for this variable.
(*bCondition)[variable] = notDefinedCondition;
returnval->AddBoundaryConditions(bregions.size()+numNewBc,bCondition);
++numNewBc;
}
}
return returnval;
}
GlobalLinSysSharedPtr Nektar::MultiRegions::DisContField1D::GetGlobalBndLinSys ( const GlobalLinSysKey mkey)

For a given key, returns the associated global linear system.

Definition at line 764 of file DisContField1D.cpp.

References ASSERTL0, ASSERTL1, Nektar::MultiRegions::eDirectFullMatrix, 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() != eDirectFullMatrix,
"Full matrix global systems are not supported for HDG "
"expansions");
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;
}
vector< bool > & Nektar::MultiRegions::DisContField1D::GetNegatedFluxNormal ( void  )

Definition at line 796 of file DisContField1D.cpp.

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

Referenced by v_AddTraceIntegral().

{
if(m_negatedFluxNormal.size() == 0)
{
Array<OneD, Array<OneD, StdRegions::StdExpansionSharedPtr> >
&elmtToTrace = m_traceMap->GetElmtToTrace();
for(int i = 0; i < GetExpSize(); ++i)
{
for(int v = 0; v < 2; ++v)
{
elmtToTrace[i][v]->as<LocalRegions::Expansion0D>();
if(vertExp->GetLeftAdjacentElementExp()->GetGeom()->GetGlobalID() != (*m_exp)[i]->GetGeom()->GetGlobalID())
{
m_negatedFluxNormal[2*i+v] = true;
}
else
{
m_negatedFluxNormal[2*i+v] = false;
}
}
}
}
}
bool Nektar::MultiRegions::DisContField1D::IsLeftAdjacentVertex ( const int  n,
const int  e 
)
private

Definition at line 253 of file DisContField1D.cpp.

References ASSERTL2, Nektar::iterator, m_boundaryVerts, Nektar::MultiRegions::ExpList::m_exp, m_periodicVerts, m_trace, and m_traceMap.

Referenced by SetUpDG().

{
m_traceMap->GetElmtToTrace()[n][e]->as<LocalRegions::Expansion0D>();
bool fwd = true;
if (traceEl->GetLeftAdjacentElementVertex () == -1 ||
traceEl->GetRightAdjacentElementVertex() == -1)
{
// Boundary edge (1 connected element). Do nothing in
// serial.
it = m_boundaryVerts.find(traceEl->GetElmtId());
// If the edge does not have a boundary condition set on
// it, then assume it is a partition edge or periodic.
if (it == m_boundaryVerts.end())
{
int traceGeomId = traceEl->GetGeom0D()->GetGlobalID();
traceGeomId);
if (pIt != m_periodicVerts.end() && !pIt->second[0].isLocal)
{
fwd = traceGeomId == min(traceGeomId,pIt->second[0].id);
}
else
{
int offset = m_trace->GetPhys_Offset(traceEl->GetElmtId());
fwd = m_traceMap->
GetTraceToUniversalMapUnique(offset) >= 0;
}
}
}
else if (traceEl->GetLeftAdjacentElementVertex () != -1 &&
traceEl->GetRightAdjacentElementVertex() != -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;
}
void Nektar::MultiRegions::DisContField1D::SetBoundaryConditionExpansion ( const SpatialDomains::MeshGraphSharedPtr graph1D,
const SpatialDomains::BoundaryConditions bcs,
const std::string  variable,
Array< OneD, MultiRegions::ExpListSharedPtr > &  bndCondExpansions,
Array< OneD, SpatialDomains::BoundaryConditionShPtr > &  bndConditions 
)
protected

Populates the list of boundary condition expansions.

Parameters
graph1DA mesh containing information about the domain and the Spectral/hp element expansion.
bcsInformation about the boundary conditions.
variableSpecifies the field.
bndCondExpansionsArray of ExpList1D objects each containing a 1D spectral/hp element expansion on a single boundary region.
bncConditionsArray of BoundaryCondition objects which contain information about the boundary conditions on the different boundary regions.

Definition at line 656 of file DisContField1D.cpp.

References ASSERTL0, Nektar::SpatialDomains::eDirichlet, Nektar::SpatialDomains::eNeumann, Nektar::SpatialDomains::eNotDefined, Nektar::SpatialDomains::ePeriodic, Nektar::SpatialDomains::eRobin, Nektar::MultiRegions::ExpList::GetBoundaryCondition(), Nektar::SpatialDomains::BoundaryConditions::GetBoundaryConditions(), Nektar::SpatialDomains::BoundaryConditions::GetBoundaryRegions(), and Nektar::iterator.

Referenced by GenerateBoundaryConditionExpansion().

{
int k;
int cnt = 0;
SpatialDomains::BoundaryRegionCollection::const_iterator it;
cnt = 0;
// list Dirichlet boundaries first
for (it = bregions.begin(); it != bregions.end(); ++it)
{
bconditions, it->first, variable);
if (locBCond->GetBoundaryConditionType() ==
{
for (bregionIt = it->second->begin();
bregionIt != it->second->end(); bregionIt++)
{
for (k = 0; k < bregionIt->second->size(); k++)
{
if ((vert = boost::dynamic_pointer_cast
(*bregionIt->second)[k])))
{
locPointExp
bndCondExpansions[cnt] = locPointExp;
bndConditions[cnt++] = locBCond;
}
else
{
ASSERTL0(false,
"dynamic cast to a vertex failed");
}
}
}
}
} // end if Dirichlet
// then, list the other (non-periodic) boundaries
for (it = bregions.begin(); it != bregions.end(); ++it)
{
locBCond = GetBoundaryCondition(bconditions, it->first, variable);
switch(locBCond->GetBoundaryConditionType())
{
case SpatialDomains::eNotDefined: // presume this will be reused as Neuman, Robin or Dirichlet later
{
for (bregionIt = it->second->begin();
bregionIt != it->second->end(); bregionIt++)
{
for (k = 0; k < bregionIt->second->size(); k++)
{
if((vert = boost::dynamic_pointer_cast
(*bregionIt->second)[k])))
{
locPointExp
bndCondExpansions[cnt] = locPointExp;
bndConditions[cnt++] = locBCond;
}
else
{
ASSERTL0(false,
"dynamic cast to a vertex failed");
}
}
}
}
// do nothing for these types
break;
default:
ASSERTL0(false,"This type of BC not implemented yet");
break;
}
}
}
void Nektar::MultiRegions::DisContField1D::SetMultiDomainBoundaryConditionExpansion ( const SpatialDomains::MeshGraphSharedPtr graph1D,
const SpatialDomains::BoundaryConditions bcs,
const std::string  variable,
Array< OneD, MultiRegions::ExpListSharedPtr > &  bndCondExpansions,
Array< OneD, SpatialDomains::BoundaryConditionShPtr > &  bndConditions,
int  subdomain 
)
protected

Populates the list of boundary condition expansions in multidomain case.

void Nektar::MultiRegions::DisContField1D::SetUpDG ( const std::string &  variable)
private

Definition at line 99 of file DisContField1D.cpp.

References ASSERTL0, ASSERTL1, Nektar::SpatialDomains::ePeriodic, Nektar::MultiRegions::_PeriodicEntity::id, IsLeftAdjacentVertex(), Nektar::MultiRegions::_PeriodicEntity::isLocal, Nektar::iterator, m_bndCondExpansions, m_bndConditions, m_boundaryVerts, Nektar::MultiRegions::ExpList::m_exp, m_globalBndMat, Nektar::MultiRegions::ExpList::m_graph, m_leftAdjacentVerts, m_periodicBwdCopy, m_periodicFwdCopy, m_periodicVerts, Nektar::MultiRegions::ExpList::m_session, m_trace, m_traceMap, Nektar::MultiRegions::NullExpListSharedPtr, Nektar::LocalRegions::Expansion0D::SetAdjacentElementExp(), and Nektar::MultiRegions::ExpList::SetUpPhysNormals().

Referenced by DisContField1D().

{
// Check for multiple calls
{
return;
}
boost::dynamic_pointer_cast<SpatialDomains::MeshGraph1D>(
*m_exp,graph1D,
m_trace = boost::dynamic_pointer_cast<ExpList>(trace);
AllocateSharedPtr(m_session, graph1D, trace, *this,
m_periodicVerts, variable);
// Scatter trace points to 1D elements. For each element, we find
// the trace point associated to each vertex. The element then
// retains a pointer to the trace space points, to ensure
// uniqueness of normals when retrieving from two adjoining
// elements which do not lie in a plane.
int ElmtPointGeom = 0;
int TracePointGeom = 0;
for (int i = 0; i < m_exp->size(); ++i)
{
exp1d = (*m_exp)[i]->as<LocalRegions::Expansion1D>();
for (int j = 0; j < exp1d->GetNverts(); ++j)
{
ElmtPointGeom = (exp1d->GetGeom1D())->GetVid(j);
for (int k = 0; k < m_trace->GetExpSize(); ++k)
{
TracePointGeom = m_trace->GetExp(k)->GetGeom()->GetVid(0);
if (TracePointGeom == ElmtPointGeom)
{
exp0d = m_trace->GetExp(k)->as<LocalRegions::Expansion0D>();
exp0d->SetAdjacentElementExp(j,exp1d);
break;
}
}
}
}
int cnt, n, e;
// Identify boundary verts
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));
}
}
else
{
ASSERTL0(false,"Periodic verts need setting up");
}
cnt += m_bndCondExpansions[n]->GetExpSize();
}
// Set up left-adjacent edge list.
m_leftAdjacentVerts.resize(2*((*m_exp).size()));
// count size of trace
for (cnt = n = 0; n < m_exp->size(); ++n)
{
for (int v = 0; v < (*m_exp)[n]->GetNverts(); ++v, ++cnt)
{
}
}
boost::unordered_map<int,pair<int,int> > perVertToExpMap;
boost::unordered_map<int,pair<int,int> >::iterator it2;
for (n = 0; n < m_exp->size(); ++n)
{
for (int v = 0; v < (*m_exp)[n]->GetNverts(); ++v)
{
(*m_exp)[n]->GetGeom()->GetVid(v));
if (it != m_periodicVerts.end())
{
perVertToExpMap[it->first] = make_pair(n,v);
}
}
}
// Set up mapping to copy Fwd of periodic bcs to Bwd of other edge.
for (n = 0; n < m_exp->size(); ++n)
{
for (int v = 0; v < (*m_exp)[n]->GetNverts(); ++v)
{
int vertGeomId = (*m_exp)[n]->GetGeom()->GetVid(v);
// Check to see if this face is periodic.
PeriodicMap::iterator it = m_periodicVerts.find(vertGeomId);
if (it != m_periodicVerts.end())
{
const PeriodicEntity &ent = it->second[0];
it2 = perVertToExpMap.find(ent.id);
if (it2 == perVertToExpMap.end())
{
if (m_session->GetComm()->GetRowComm()->GetSize() > 1 &&
!ent.isLocal)
{
continue;
}
else
{
ASSERTL1(false, "Periodic vert not found!");
}
}
int offset = m_trace->GetPhys_Offset((m_traceMap->GetElmtToTrace())
[n][v]->GetElmtId());
m_periodicFwdCopy.push_back(offset);
int offset2 = m_trace->GetPhys_Offset((m_traceMap->GetElmtToTrace())
[it2->second.first]
[it2->second.second]->GetElmtId());
m_periodicBwdCopy.push_back(offset2);
}
}
}
}
void Nektar::MultiRegions::DisContField1D::v_AddTraceIntegral ( const Array< OneD, const NekDouble > &  Fn,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1028 of file DisContField1D.cpp.

References Nektar::LibUtilities::BasisManager(), Nektar::LibUtilities::eGauss_Lagrange, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::MultiRegions::ExpList::GetCoeff_Offset(), Nektar::MultiRegions::ExpList::GetExpSize(), GetNegatedFluxNormal(), Nektar::MultiRegions::ExpList::GetTrace(), and m_traceMap.

{
int n,offset, t_offset;
Array<OneD, Array<OneD, StdRegions::StdExpansionSharedPtr> >
&elmtToTrace = m_traceMap->GetElmtToTrace();
// Basis shared pointer
vector<bool> negatedFluxNormal = GetNegatedFluxNormal();
for (n = 0; n < GetExpSize(); ++n)
{
// Basis definition on each element
Basis = (*m_exp)[n]->GetBasis(0);
// Number of coefficients on each element
int e_ncoeffs = (*m_exp)[n]->GetNcoeffs();
offset = GetCoeff_Offset(n);
// Implementation for every points except Gauss points
if (Basis->GetBasisType() != LibUtilities::eGauss_Lagrange)
{
t_offset = GetTrace()->GetCoeff_Offset(elmtToTrace[n][0]->GetElmtId());
if(negatedFluxNormal[2*n])
{
outarray[offset] -= Fn[t_offset];
}
else
{
outarray[offset] += Fn[t_offset];
}
t_offset = GetTrace()->GetCoeff_Offset(elmtToTrace[n][1]->GetElmtId());
if(negatedFluxNormal[2*n+1])
{
outarray[offset+(*m_exp)[n]->GetVertexMap(1)] -= Fn[t_offset];
}
else
{
outarray[offset+(*m_exp)[n]->GetVertexMap(1)] += Fn[t_offset];
}
}
else
{
#if 0
BS_k(LibUtilities::eGauss_Lagrange,e_ncoeffs,BS_p);
Array<OneD, NekDouble> coords(3, 0.0);
int j;
for(p = 0; p < 2; ++p)
{
NekDouble vertnorm = 0.0;
for (int i=0; i<((*m_exp)[n]->
GetVertexNormal(p)).num_elements(); i++)
{
vertnorm += ((*m_exp)[n]->GetVertexNormal(p))[i][0];
coords[0] = vertnorm ;
}
t_offset = GetTrace()->GetPhys_Offset(n+p);
if (vertnorm >= 0.0)
{
m_Ixm = BASE->GetI(coords);
for (j = 0; j < e_ncoeffs; j++)
{
outarray[offset + j] +=
(m_Ixm->GetPtr())[j] * Fn[t_offset];
}
}
if (vertnorm < 0.0)
{
m_Ixm = BASE->GetI(coords);
for (j = 0; j < e_ncoeffs; j++)
{
outarray[offset + j] -=
(m_Ixm->GetPtr())[j] * Fn[t_offset];
}
}
}
#else
int j;
static DNekMatSharedPtr m_Ixm, m_Ixp;
static int sav_ncoeffs = 0;
if(!m_Ixm.get() || e_ncoeffs != sav_ncoeffs)
{
BS_k(LibUtilities::eGauss_Lagrange,e_ncoeffs,BS_p);
Array<OneD, NekDouble> coords(1, 0.0);
coords[0] = -1.0;
m_Ixm = BASE->GetI(coords);
coords[0] = 1.0;
m_Ixp = BASE->GetI(coords);
sav_ncoeffs = e_ncoeffs;
}
t_offset = GetTrace()->GetCoeff_Offset(elmtToTrace[n][0]->GetElmtId());
if(negatedFluxNormal[2*n])
{
for (j = 0; j < e_ncoeffs; j++)
{
outarray[offset + j] -=
(m_Ixm->GetPtr())[j] * Fn[t_offset];
}
}
else
{
for (j = 0; j < e_ncoeffs; j++)
{
outarray[offset + j] +=
(m_Ixm->GetPtr())[j] * Fn[t_offset];
}
}
t_offset = GetTrace()->GetCoeff_Offset(elmtToTrace[n][1]->GetElmtId());
if (negatedFluxNormal[2*n+1])
{
for (j = 0; j < e_ncoeffs; j++)
{
outarray[offset + j] -=
(m_Ixp->GetPtr())[j] * Fn[t_offset];
}
}
else
{
for (j = 0; j < e_ncoeffs; j++)
{
outarray[offset + j] +=
(m_Ixp->GetPtr())[j] * Fn[t_offset];
}
}
#endif
}
}
}
void Nektar::MultiRegions::DisContField1D::v_EvaluateBoundaryConditions ( const NekDouble  time = 0.0,
const std::string  varName = "",
const NekDouble  x2_in = NekConstants::kNekUnsetDouble,
const NekDouble  x3_in = NekConstants::kNekUnsetDouble 
)
protectedvirtual

Evaluate all boundary conditions at a given time..

Based on the expression $g(x,t)$ for the boundary conditions, this function evaluates the boundary conditions for all boundaries at time-level t.

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

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1322 of file DisContField1D.cpp.

References ASSERTL0, Nektar::SpatialDomains::eDirichlet, Nektar::SpatialDomains::eNeumann, Nektar::SpatialDomains::eNotDefined, Nektar::SpatialDomains::eQinflow, Nektar::SpatialDomains::eRCRterminal, Nektar::SpatialDomains::eRobin, Nektar::SpatialDomains::eTimeDependent, Nektar::MultiRegions::ExpList::GetCoeff(), Nektar::NekConstants::kNekUnsetDouble, m_bndCondExpansions, and m_bndConditions.

{
int i;
Array<OneD, NekDouble> x0(1);
Array<OneD, NekDouble> x1(1);
Array<OneD, NekDouble> x2(1);
for (i = 0; i < m_bndCondExpansions.num_elements(); ++i)
{
if (time == 0.0 || m_bndConditions[i]->GetUserDefined() ==
m_bndConditions[i]->GetUserDefined() ==
m_bndConditions[i]->GetUserDefined() ==
{
m_bndCondExpansions[i]->GetCoords(x0, x1, x2);
if (x2_in != NekConstants::kNekUnsetDouble && x3_in !=
{
x1[0] = x2_in;
x2[0] = x3_in;
}
if (m_bndConditions[i]->GetBoundaryConditionType() ==
{
m_bndCondExpansions[i]->SetCoeff(0,
(boost::static_pointer_cast<SpatialDomains
::DirichletBoundaryCondition>(m_bndConditions[i])
->m_dirichletCondition).Evaluate(x0[0],x1[0],x2[0],time));
}
else if (m_bndConditions[i]->GetBoundaryConditionType()
{
m_bndCondExpansions[i]->SetCoeff(0,
(boost::static_pointer_cast<SpatialDomains
::NeumannBoundaryCondition>(m_bndConditions[i])
->m_neumannCondition).Evaluate(x0[0],x1[0],x2[0],time));
}
else if (m_bndConditions[i]->GetBoundaryConditionType()
{
m_bndCondExpansions[i]->SetCoeff(0,
(boost::static_pointer_cast<SpatialDomains
::RobinBoundaryCondition>(m_bndConditions[i])
->m_robinFunction).Evaluate(x0[0],x1[0],x2[0],time));
m_bndCondExpansions[i]->SetPhys(0,
(boost::static_pointer_cast<SpatialDomains
::RobinBoundaryCondition>(m_bndConditions[i])
->m_robinPrimitiveCoeff).Evaluate(x0[0],x1[0],x2[0],time));
}
else if (m_bndConditions[i]->GetBoundaryConditionType()
{
}
else
{
ASSERTL0(false, "This type of BC not implemented yet");
}
}
}
}
void Nektar::MultiRegions::DisContField1D::v_ExtractTracePhys ( Array< OneD, NekDouble > &  outarray)
protectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 983 of file DisContField1D.cpp.

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

{
ASSERTL1(m_physState == true,"local physical space is not true ");
}
void Nektar::MultiRegions::DisContField1D::v_ExtractTracePhys ( const Array< OneD, const NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray 
)
protectedvirtual

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

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

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

This will not work for non-boundary expansions

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1003 of file DisContField1D.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,p,offset,phys_offset;
ASSERTL1(outarray.num_elements() >= m_trace->GetExpSize(),
"input array is of insufficient length");
for (n = 0; n < nexp; ++n)
{
phys_offset = GetPhys_Offset(n);
for (p = 0; p < (*m_exp)[n]->GetNverts(); ++p)
{
offset = m_trace->GetPhys_Offset(
(m_traceMap->GetElmtToTrace())[n][p]->GetElmtId());
(*m_exp)[n]->GetVertexPhysVals(p, inarray + phys_offset,
outarray[offset]);
}
}
}
virtual const Array<OneD,const MultiRegions::ExpListSharedPtr>& Nektar::MultiRegions::DisContField1D::v_GetBndCondExpansions ( void  )
inlineprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 218 of file DisContField1D.h.

References m_bndCondExpansions.

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

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 224 of file DisContField1D.h.

References m_bndConditions.

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

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1397 of file DisContField1D.cpp.

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

{
map<int, int> VertGID;
int i,n,id;
int bid,cnt,Vid;
int nbcs = m_bndConditions.num_elements();
// make sure arrays are of sufficient length
if (ElmtID.num_elements() != nbcs)
{
ElmtID = Array<OneD, int>(nbcs,-1);
}
else
{
fill(ElmtID.get(), ElmtID.get()+nbcs, -1);
}
if (VertID.num_elements() != nbcs)
{
VertID = Array<OneD, int>(nbcs);
}
// setup map of all global ids along boundary
for (cnt = n = 0; n < m_bndCondExpansions.num_elements(); ++n)
{
Vid = m_bndCondExpansions[n]->GetExp(0)->GetGeom()->GetVertex(0)->GetVid();
VertGID[Vid] = cnt++;
}
// Loop over elements and find verts that match;
for(cnt = n = 0; n < GetExpSize(); ++n)
{
exp = (*m_exp)[n];
for(i = 0; i < exp->GetNverts(); ++i)
{
id = exp->GetGeom()->GetVid(i);
if (VertGID.count(id) > 0)
{
bid = VertGID.find(id)->second;
ASSERTL1(ElmtID[bid] == -1,"Edge already set");
ElmtID[bid] = n;
VertID[bid] = i;
cnt ++;
}
}
}
ASSERTL1(cnt == nbcs,"Failed to visit all boundary condtiions");
}
void Nektar::MultiRegions::DisContField1D::v_GetFwdBwdTracePhys ( Array< OneD, NekDouble > &  Fwd,
Array< OneD, NekDouble > &  Bwd 
)
protectedvirtual

Generate the forward or backward state for each trace point.

Parameters
FwdForward state.
BwdBackward state.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 835 of file DisContField1D.cpp.

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

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

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

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

If the element is left-adjacent we extract the vertex 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 1D data from which we wish to extract the backward and forward orientated trace/edge arrays.
FwdThe resulting forwards space.
BwdThe resulting backwards space.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 868 of file DisContField1D.cpp.

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

{
// Expansion casts
// Counter variables
int n, v;
// Number of elements
int nElements = GetExpSize();
// Number of solution points of each element
int nLocalSolutionPts;
// Initial index of each element
int phys_offset;
Array<OneD, Array<OneD, StdRegions::StdExpansionSharedPtr> >
&elmtToTrace = m_traceMap->GetElmtToTrace();
// Basis shared pointer
// Set forward and backard state to zero
Vmath::Zero(Fwd.num_elements(), Fwd, 1);
Vmath::Zero(Bwd.num_elements(), Bwd, 1);
int cnt;
// Loop on the elements
for (cnt = n = 0; n < nElements; ++n)
{
exp1D = (*m_exp)[n]->as<LocalRegions::Expansion1D>();
// Set the offset of each element
phys_offset = GetPhys_Offset(n);
// Set the number of solution points of each element
nLocalSolutionPts = (*m_exp)[n]->GetNumPoints(0);
Basis = (*m_exp)[n]->GetBasis(0);
for(v = 0; v < 2; ++v, ++cnt)
{
int offset = m_trace->GetPhys_Offset(elmtToTrace[n][v]->GetElmtId());
{
(*m_exp)[n]->GetVertexPhysVals(v, field + phys_offset,
Fwd[offset]);
}
else
{
(*m_exp)[n]->GetVertexPhysVals(v, field + phys_offset,
Bwd[offset]);
}
}
}
// Fill boundary conditions into missing elements.
int id = 0;
for(cnt = n = 0; n < m_bndCondExpansions.num_elements(); ++n)
{
if (m_bndConditions[n]->GetBoundaryConditionType() ==
{
id = m_trace->GetPhys_Offset(m_traceMap->GetBndCondTraceToGlobalTraceMap(cnt));
Bwd[id] = m_bndCondExpansions[n]->GetPhys()[0]; //this is not getting the correct value?
cnt++;
}
else if (m_bndConditions[n]->GetBoundaryConditionType() ==
m_bndConditions[n]->GetBoundaryConditionType() ==
{
"Method not set up for non-zero Neumann "
"boundary condition");
id = m_trace->GetPhys_Offset(m_traceMap->GetBndCondTraceToGlobalTraceMap(cnt));
Bwd[id] = Fwd[id];
cnt++;
}
else if (m_bndConditions[n]->GetBoundaryConditionType() ==
{
// Do nothing
}
else if (m_bndConditions[n]->GetBoundaryConditionType() !=
{
ASSERTL0(false,
"Method not set up for this boundary condition.");
}
}
// Copy any periodic boundary conditions.
for (n = 0; n < m_periodicFwdCopy.size(); ++n)
{
}
// Do parallel exchange for forwards/backwards spaces.
m_traceMap->UniversalTraceAssemble(Fwd);
m_traceMap->UniversalTraceAssemble(Bwd);
}
map< int, RobinBCInfoSharedPtr > Nektar::MultiRegions::DisContField1D::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 1464 of file DisContField1D.cpp.

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

{
int i;
map<int, RobinBCInfoSharedPtr> returnval;
Array<OneD, int> ElmtID,VertID;
GetBoundaryToElmtMap(ElmtID,VertID);
for (i = 0; i < m_bndCondExpansions.num_elements(); ++i)
{
if (m_bndConditions[i]->GetBoundaryConditionType() ==
{
int elmtid;
Array<OneD, NekDouble> Array_tmp;
locExpList = m_bndCondExpansions[i];
VertID[i],Array_tmp = locExpList->GetPhys());
elmtid = ElmtID[i];
// make link list if necessary (not likely in
// 1D but needed in 2D & 3D)
if(returnval.count(elmtid) != 0)
{
rInfo->next = returnval.find(elmtid)->second;
}
returnval[elmtid] = rInfo;
}
}
return returnval;
}
virtual ExpListSharedPtr& Nektar::MultiRegions::DisContField1D::v_GetTrace ( )
inlineprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 164 of file DisContField1D.h.

References m_trace.

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

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 169 of file DisContField1D.h.

References m_traceMap.

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

Solve the Helmholtz equation.

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 1194 of file DisContField1D.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::IProductWRTBase(), m_bndCondExpansions, m_bndConditions, Nektar::MultiRegions::ExpList::m_coeff_offset, Nektar::MultiRegions::ExpList::m_ncoeffs, m_traceMap, Vmath::Neg(), Nektar::MultiRegions::NullAssemblyMapSharedPtr, Nektar::Transpose(), and Vmath::Zero().

{
int i,n,cnt,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 Boundary System
//----------------------------------
int GloBndDofs = m_traceMap->GetNumGlobalBndCoeffs();
int NumDirBCs = m_traceMap->GetNumLocalDirBndCoeffs();
int e_ncoeffs,id;
GlobalMatrixKey HDGLamToUKey(
factors,
varcoeff);
const DNekScalBlkMatSharedPtr &HDGLamToU =
GetBlockMatrix(HDGLamToUKey);
// Retrieve global trace space storage, \Lambda, from trace expansion
Array<OneD,NekDouble> BndSol = Array<OneD,NekDouble>
(m_traceMap->GetNumLocalBndCoeffs());
Array<OneD,NekDouble> BndRhs(GloBndDofs,0.0);
// Zero trace space
Vmath::Zero(GloBndDofs,BndSol,1);
int LocBndCoeffs = m_traceMap->GetNumLocalBndCoeffs();
Array<OneD, NekDouble> loc_lambda(LocBndCoeffs);
DNekVec LocLambda(LocBndCoeffs,loc_lambda,eWrapper);
//----------------------------------
// Evaluate Trace Forcing
//----------------------------------
// Determing <u_lam,f> terms using HDGLamToU matrix
for (cnt = n = 0; n < nexp; ++n)
{
nbndry = (*m_exp)[n]->NumDGBndryCoeffs();
e_ncoeffs = (*m_exp)[n]->GetNcoeffs();
e_f = f+m_coeff_offset[n];
e_l = loc_lambda + cnt;
// use outarray as tmp space
DNekVec Floc (nbndry, e_l, eWrapper);
DNekVec ElmtFce (e_ncoeffs, e_f, eWrapper);
Floc = Transpose(*(HDGLamToU->GetBlock(n,n)))*ElmtFce;
cnt += nbndry;
}
// Assemble into global operator
m_traceMap->AssembleBnd(loc_lambda,BndRhs);
cnt = 0;
// Copy Dirichlet boundary conditions into trace space
for (i = 0; i < m_bndCondExpansions.num_elements(); ++i)
{
if (m_bndConditions[i]->GetBoundaryConditionType() ==
{
id = m_traceMap->GetBndCondCoeffsToGlobalCoeffsMap(i);
BndSol[id] = m_bndCondExpansions[i]->GetCoeff(0);
}
else
{
id = m_traceMap->GetBndCondCoeffsToGlobalCoeffsMap(i);
BndRhs[id] += m_bndCondExpansions[i]->GetCoeff(0);
}
}
//----------------------------------
// Solve trace problem
//----------------------------------
if (GloBndDofs - NumDirBCs > 0)
{
GlobalLinSysKey key(StdRegions::eHybridDGHelmBndLam,
m_traceMap,factors);
LinSys->Solve(BndRhs,BndSol,m_traceMap);
}
//----------------------------------
// Internal element solves
//----------------------------------
GlobalMatrixKey invHDGhelmkey(StdRegions::eInvHybridDGHelmholtz,
factors);
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::DisContField1D::v_UpdateBndCondExpansion ( int  i)
inlineprotectedvirtual

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 230 of file DisContField1D.h.

References m_bndCondExpansions.

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

Reimplemented from Nektar::MultiRegions::ExpList.

Definition at line 236 of file DisContField1D.h.

References m_bndConditions.

{
}

Member Data Documentation

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

Discretised boundary conditions.

It is an array of size equal to the number of boundary points and consists of entries of the type LocalRegions::PointExp. Every entry corresponds to a point on a single boundary region.

Definition at line 110 of file DisContField1D.h.

Referenced by Nektar::MultiRegions::ContField1D::ContField1D(), GenerateBoundaryConditionExpansion(), Nektar::MultiRegions::ContField1D::GetBndCondExpansions(), Nektar::MultiRegions::ContField1D::GlobalSolve(), SetUpDG(), v_EvaluateBoundaryConditions(), v_GetBndCondExpansions(), v_GetBoundaryToElmtMap(), v_GetFwdBwdTracePhys(), v_GetRobinBCInfo(), Nektar::MultiRegions::ContField1D::v_HelmSolve(), v_HelmSolve(), Nektar::MultiRegions::ContField1D::v_ImposeDirichletConditions(), and v_UpdateBndCondExpansion().

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

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

Definition at line 114 of file DisContField1D.h.

Referenced by Nektar::MultiRegions::ContField1D::ContField1D(), GenerateBoundaryConditionExpansion(), Nektar::MultiRegions::ContField1D::GetBndConditions(), Nektar::MultiRegions::ContField1D::GlobalSolve(), SetUpDG(), v_EvaluateBoundaryConditions(), v_GetBndConditions(), v_GetBoundaryToElmtMap(), v_GetFwdBwdTracePhys(), v_GetRobinBCInfo(), Nektar::MultiRegions::ContField1D::v_HelmSolve(), v_HelmSolve(), Nektar::MultiRegions::ContField1D::v_ImposeDirichletConditions(), and v_UpdateBndConditions().

std::set<int> Nektar::MultiRegions::DisContField1D::m_boundaryVerts
protected

A set storing the global IDs of any boundary edges.

Definition at line 128 of file DisContField1D.h.

Referenced by IsLeftAdjacentVertex(), and SetUpDG().

GlobalLinSysMapShPtr Nektar::MultiRegions::DisContField1D::m_globalBndMat
protected

Global boundary matrix.

Definition at line 117 of file DisContField1D.h.

Referenced by GetGlobalBndLinSys(), and SetUpDG().

vector<bool> Nektar::MultiRegions::DisContField1D::m_leftAdjacentVerts
protected

Definition at line 150 of file DisContField1D.h.

Referenced by SetUpDG(), and v_GetFwdBwdTracePhys().

vector<bool> Nektar::MultiRegions::DisContField1D::m_negatedFluxNormal
private

Definition at line 265 of file DisContField1D.h.

Referenced by GetNegatedFluxNormal().

int Nektar::MultiRegions::DisContField1D::m_numDirBndCondExpansions
protected

The number of boundary segments on which Dirichlet boundary conditions are imposed.

Definition at line 102 of file DisContField1D.h.

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

Definition at line 143 of file DisContField1D.h.

Referenced by SetUpDG(), and v_GetFwdBwdTracePhys().

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

Referenced by SetUpDG(), and v_GetFwdBwdTracePhys().

PeriodicMap Nektar::MultiRegions::DisContField1D::m_periodicVerts
protected

A map which identifies groups of periodic vertices.

Definition at line 134 of file DisContField1D.h.

Referenced by Nektar::MultiRegions::ContField1D::ContField1D(), FindPeriodicVertices(), IsLeftAdjacentVertex(), and SetUpDG().

ExpListSharedPtr Nektar::MultiRegions::DisContField1D::m_trace
protected

Trace space storage for points between elements.

Definition at line 120 of file DisContField1D.h.

Referenced by IsLeftAdjacentVertex(), SetUpDG(), v_ExtractTracePhys(), v_GetFwdBwdTracePhys(), and v_GetTrace().

AssemblyMapDGSharedPtr Nektar::MultiRegions::DisContField1D::m_traceMap
protected

Local to global DG mapping for trace space.

Definition at line 123 of file DisContField1D.h.

Referenced by GetGlobalBndLinSys(), GetNegatedFluxNormal(), IsLeftAdjacentVertex(), SetUpDG(), v_AddTraceIntegral(), v_ExtractTracePhys(), v_GetFwdBwdTracePhys(), v_GetTraceMap(), and v_HelmSolve().