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

Constructs mappings for the C0 scalar continuous Galerkin formulation. More...

#include <AssemblyMapCG.h>

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

Public Member Functions

 AssemblyMapCG (const LibUtilities::SessionReaderSharedPtr &pSession, const std::string variable="DefaultVar")
 Default constructor.
 AssemblyMapCG (const LibUtilities::SessionReaderSharedPtr &pSession, const int numLocalCoeffs, const ExpList &locExp)
 General constructor for expansions of all dimensions without boundary conditions.
virtual ~AssemblyMapCG ()
 Destructor.
map< int, vector< ExtraDirDof > > & GetExtraDirDofs ()
- Public Member Functions inherited from Nektar::MultiRegions::AssemblyMap
 AssemblyMap ()
 Default constructor.
 AssemblyMap (const LibUtilities::SessionReaderSharedPtr &pSession, const std::string variable="DefaultVar")
 Constructor with a communicator.
 AssemblyMap (AssemblyMap *oldLevelMap, const BottomUpSubStructuredGraphSharedPtr &multiLevelGraph)
 Constructor for next level in multi-level static condensation.
virtual ~AssemblyMap ()
 Destructor.
LibUtilities::CommSharedPtr GetComm ()
 Retrieves the communicator.
size_t GetHash () const
 Retrieves the hash of this map.
int GetLocalToGlobalMap (const int i) const
int GetGlobalToUniversalMap (const int i) const
int GetGlobalToUniversalMapUnique (const int i) const
const Array< OneD, const int > & GetLocalToGlobalMap ()
const Array< OneD, const int > & GetGlobalToUniversalMap ()
const Array< OneD, const int > & GetGlobalToUniversalMapUnique ()
NekDouble GetLocalToGlobalSign (const int i) const
const Array< OneD, NekDouble > & GetLocalToGlobalSign () const
void LocalToGlobal (const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
void LocalToGlobal (const NekVector< NekDouble > &loc, NekVector< NekDouble > &global) const
void GlobalToLocal (const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const
void GlobalToLocal (const NekVector< NekDouble > &global, NekVector< NekDouble > &loc) const
void Assemble (const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
void Assemble (const NekVector< NekDouble > &loc, NekVector< NekDouble > &global) const
void UniversalAssemble (Array< OneD, NekDouble > &pGlobal) const
void UniversalAssemble (NekVector< NekDouble > &pGlobal) const
void UniversalAssemble (Array< OneD, NekDouble > &pGlobal, int offset) const
int GetLocalToGlobalBndMap (const int i) const
 Retrieve the global index of a given local boundary mode.
const Array< OneD, const int > & GetLocalToGlobalBndMap ()
 Retrieve the global indices of the local boundary modes.
const Array< OneD, const int > & GetGlobalToUniversalBndMap ()
const Array< OneD, const int > & GetGlobalToUniversalBndMapUnique ()
bool GetSignChange ()
 Returns true if using a modal expansion requiring a change of sign of some modes.
NekDouble GetLocalToGlobalBndSign (const int i) const
 Retrieve the sign change of a given local boundary mode.
Array< OneD, const NekDoubleGetLocalToGlobalBndSign () const
 Retrieve the sign change for all local boundary modes.
int GetBndCondCoeffsToGlobalCoeffsMap (const int i)
 Retrieves the global index corresponding to a boundary expansion mode.
const Array< OneD, const int > & GetBndCondCoeffsToGlobalCoeffsMap ()
 Retrieves the global indices corresponding to the boundary expansion modes.
NekDouble GetBndCondCoeffsToGlobalCoeffsSign (const int i)
 Returns the modal sign associated with a given boundary expansion mode.
int GetBndCondTraceToGlobalTraceMap (const int i)
 Returns the global index of the boundary trace giving the index on the boundary expansion.
const Array< OneD, const int > & GetBndCondTraceToGlobalTraceMap ()
int GetNumGlobalDirBndCoeffs () const
 Returns the number of global Dirichlet boundary coefficients.
int GetNumLocalDirBndCoeffs () const
 Returns the number of local Dirichlet boundary coefficients.
int GetNumGlobalBndCoeffs () const
 Returns the total number of global boundary coefficients.
int GetNumLocalBndCoeffs () const
 Returns the total number of local boundary coefficients.
int GetNumLocalCoeffs () const
 Returns the total number of local coefficients.
int GetNumGlobalCoeffs () const
 Returns the total number of global coefficients.
bool GetSingularSystem () const
 Retrieves if the system is singular (true) or not (false)
void GlobalToLocalBnd (const NekVector< NekDouble > &global, NekVector< NekDouble > &loc, int offset) const
void GlobalToLocalBnd (const NekVector< NekDouble > &global, NekVector< NekDouble > &loc) const
void GlobalToLocalBnd (const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc, int offset) const
void GlobalToLocalBnd (const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const
void LocalBndToGlobal (const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const
void LocalBndToGlobal (const NekVector< NekDouble > &loc, NekVector< NekDouble > &global) const
void LocalBndToGlobal (const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, int offset) const
void LocalBndToGlobal (const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
void AssembleBnd (const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const
void AssembleBnd (const NekVector< NekDouble > &loc, NekVector< NekDouble > &global) const
void AssembleBnd (const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, int offset) const
void AssembleBnd (const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
void UniversalAssembleBnd (Array< OneD, NekDouble > &pGlobal) const
void UniversalAssembleBnd (NekVector< NekDouble > &pGlobal) const
void UniversalAssembleBnd (Array< OneD, NekDouble > &pGlobal, int offset) const
int GetFullSystemBandWidth () const
int GetNumNonDirVertexModes () const
int GetNumNonDirEdgeModes () const
int GetNumNonDirFaceModes () const
int GetNumDirEdges () const
int GetNumDirFaces () const
int GetNumNonDirEdges () const
int GetNumNonDirFaces () const
void PrintStats (std::ostream &out, std::string variable) const
const Array< OneD, const int > & GetExtraDirEdges ()
boost::shared_ptr< AssemblyMapLinearSpaceMap (const ExpList &locexp, GlobalSysSolnType solnType)
int GetBndSystemBandWidth () const
 Returns the bandwidth of the boundary system.
int GetStaticCondLevel () const
 Returns the level of static condensation for this map.
int GetNumPatches () const
 Returns the number of patches in this static condensation level.
const Array< OneD, const
unsigned int > & 
GetNumLocalBndCoeffsPerPatch ()
 Returns the number of local boundary coefficients in each patch.
const Array< OneD, const
unsigned int > & 
GetNumLocalIntCoeffsPerPatch ()
 Returns the number of local interior coefficients in each patch.
const AssemblyMapSharedPtr GetNextLevelLocalToGlobalMap () const
 Returns the local to global mapping for the next level in the multi-level static condensation.
void SetNextLevelLocalToGlobalMap (AssemblyMapSharedPtr pNextLevelLocalToGlobalMap)
const PatchMapSharedPtrGetPatchMapFromPrevLevel (void) const
 Returns the patch map from the previous level of the multi-level static condensation.
bool AtLastLevel () const
 Returns true if this is the last level in the multi-level static condensation.
GlobalSysSolnType GetGlobalSysSolnType () const
 Returns the method of solving global systems.
PreconditionerType GetPreconType () const
NekDouble GetIterativeTolerance () const
int GetSuccessiveRHS () const
int GetLowestStaticCondLevel () const

Protected Member Functions

void SetUpUniversalC0ContMap (const ExpList &locExp, const PeriodicMap &perVerts=NullPeriodicMap, const PeriodicMap &perEdges=NullPeriodicMap, const PeriodicMap &perFaces=NullPeriodicMap)
void CalculateFullSystemBandWidth ()
 Calculate the bandwith of the full matrix system.
virtual int v_GetLocalToGlobalMap (const int i) const
virtual int v_GetGlobalToUniversalMap (const int i) const
virtual int v_GetGlobalToUniversalMapUnique (const int i) const
virtual const Array< OneD,
const int > & 
v_GetLocalToGlobalMap ()
virtual const Array< OneD,
const int > & 
v_GetGlobalToUniversalMap ()
virtual const Array< OneD,
const int > & 
v_GetGlobalToUniversalMapUnique ()
virtual NekDouble v_GetLocalToGlobalSign (const int i) const
virtual const Array< OneD,
NekDouble > & 
v_GetLocalToGlobalSign () const
virtual void v_LocalToGlobal (const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
virtual void v_LocalToGlobal (const NekVector< NekDouble > &loc, NekVector< NekDouble > &global) const
virtual void v_GlobalToLocal (const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const
virtual void v_GlobalToLocal (const NekVector< NekDouble > &global, NekVector< NekDouble > &loc) const
virtual void v_Assemble (const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
virtual void v_Assemble (const NekVector< NekDouble > &loc, NekVector< NekDouble > &global) const
virtual void v_UniversalAssemble (Array< OneD, NekDouble > &pGlobal) const
virtual void v_UniversalAssemble (NekVector< NekDouble > &pGlobal) const
virtual void v_UniversalAssemble (Array< OneD, NekDouble > &pGlobal, int offset) const
virtual int v_GetFullSystemBandWidth () const
virtual int v_GetNumNonDirVertexModes () const
virtual int v_GetNumNonDirEdgeModes () const
virtual int v_GetNumNonDirFaceModes () const
virtual int v_GetNumDirEdges () const
virtual int v_GetNumDirFaces () const
virtual int v_GetNumNonDirEdges () const
virtual int v_GetNumNonDirFaces () const
virtual const Array< OneD,
const int > & 
v_GetExtraDirEdges ()
virtual AssemblyMapSharedPtr v_LinearSpaceMap (const ExpList &locexp, GlobalSysSolnType solnType)
 Construct an AssemblyMapCG object which corresponds to the linear space of the current object.
- Protected Member Functions inherited from Nektar::MultiRegions::AssemblyMap
void CalculateBndSystemBandWidth ()
 Calculates the bandwidth of the boundary system.
void GlobalToLocalBndWithoutSign (const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc)

Protected Attributes

Array< OneD, int > m_localToGlobalMap
 Integer map of local coeffs to global space.
Array< OneD, NekDoublem_localToGlobalSign
 Integer sign of local coeffs to global space.
int m_fullSystemBandWidth
 Bandwith of the full matrix system (no static condensation).
Array< OneD, int > m_globalToUniversalMap
 Integer map of process coeffs to universal space.
Array< OneD, int > m_globalToUniversalMapUnique
 Integer map of unique process coeffs to universal space (signed)
int m_numNonDirVertexModes
 Number of non Dirichlet vertex modes.
int m_numNonDirEdgeModes
 Number of non Dirichlet edge modes.
int m_numNonDirFaceModes
 Number of non Dirichlet face modes.
int m_numDirEdges
 Number of Dirichlet edges.
int m_numDirFaces
 Number of Dirichlet faces.
int m_numNonDirEdges
 Number of Dirichlet edges.
int m_numNonDirFaces
 Number of Dirichlet faces.
Array< OneD, int > m_extraDirEdges
 Extra dirichlet edges in parallel.
int m_maxStaticCondLevel
 Maximum static condensation level.
map< int, vector< ExtraDirDof > > m_extraDirDofs
 Map indicating degrees of freedom which are Dirichlet but whose value is stored on another processor.
- Protected Attributes inherited from Nektar::MultiRegions::AssemblyMap
LibUtilities::SessionReaderSharedPtr m_session
 Session object.
LibUtilities::CommSharedPtr m_comm
 Communicator.
size_t m_hash
 Hash for map.
int m_numLocalBndCoeffs
 Number of local boundary coefficients.
int m_numGlobalBndCoeffs
 Total number of global boundary coefficients.
int m_numLocalDirBndCoeffs
 Number of Local Dirichlet Boundary Coefficients.
int m_numGlobalDirBndCoeffs
 Number of Global Dirichlet Boundary Coefficients.
bool m_systemSingular
 Flag indicating if the system is singular or not.
int m_numLocalCoeffs
 Total number of local coefficients.
int m_numGlobalCoeffs
 Total number of global coefficients.
bool m_signChange
 Flag indicating if modes require sign reversal.
Array< OneD, int > m_localToGlobalBndMap
 Integer map of local boundary coeffs to global space.
Array< OneD, NekDoublem_localToGlobalBndSign
 Integer sign of local boundary coeffs to global space.
Array< OneD, int > m_bndCondCoeffsToGlobalCoeffsMap
 Integer map of bnd cond coeffs to global coefficients.
Array< OneD, NekDoublem_bndCondCoeffsToGlobalCoeffsSign
 Integer map of bnd cond coeffs to global coefficients.
Array< OneD, int > m_bndCondTraceToGlobalTraceMap
 Integer map of bnd cond trace number to global trace number.
Array< OneD, int > m_globalToUniversalBndMap
 Integer map of process coeffs to universal space.
Array< OneD, int > m_globalToUniversalBndMapUnique
 Integer map of unique process coeffs to universal space (signed)
GlobalSysSolnType m_solnType
 The solution type of the global system.
int m_bndSystemBandWidth
 The bandwith of the global bnd system.
PreconditionerType m_preconType
 Type type of preconditioner to use in iterative solver.
NekDouble m_iterativeTolerance
 Tolerance for iterative solver.
int m_successiveRHS
 sucessive RHS for iterative solver
Gs::gs_datam_gsh
Gs::gs_datam_bndGsh
int m_staticCondLevel
 The level of recursion in the case of multi-level static condensation.
int m_numPatches
 The number of patches (~elements) in the current level.
Array< OneD, unsigned int > m_numLocalBndCoeffsPerPatch
 The number of bnd dofs per patch.
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
 The number of int dofs per patch.
AssemblyMapSharedPtr m_nextLevelLocalToGlobalMap
 Map from the patches of the previous level to the patches of the current level.
int m_lowestStaticCondLevel
 Lowest static condensation level.

Detailed Description

Constructs mappings for the C0 scalar continuous Galerkin formulation.

Mappings are created for three possible global solution types:

These mappings are used by GlobalLinSys to generate the global system.

Definition at line 63 of file AssemblyMapCG.h.

Constructor & Destructor Documentation

Nektar::MultiRegions::AssemblyMapCG::AssemblyMapCG ( const LibUtilities::SessionReaderSharedPtr pSession,
const std::string  variable = "DefaultVar" 
)

Default constructor.

Definition at line 69 of file AssemblyMapCG.cpp.

References m_maxStaticCondLevel.

:
AssemblyMap(pSession,variable)
{
pSession->LoadParameter("MaxStaticCondLevel",m_maxStaticCondLevel,100);
}
Nektar::MultiRegions::AssemblyMapCG::AssemblyMapCG ( const LibUtilities::SessionReaderSharedPtr pSession,
const int  numLocalCoeffs,
const ExpList locExp 
)

General constructor for expansions of all dimensions without boundary conditions.

Definition at line 81 of file AssemblyMapCG.cpp.

References ASSERTL0.

:
AssemblyMap(pSession)
{
ASSERTL0(false,"AssemblyMapCG: you need to instantiate dimension-specific derived class.");
}
Nektar::MultiRegions::AssemblyMapCG::~AssemblyMapCG ( )
virtual

Destructor.

Definition at line 95 of file AssemblyMapCG.cpp.

{
}

Member Function Documentation

void Nektar::MultiRegions::AssemblyMapCG::CalculateFullSystemBandWidth ( )
protected

Calculate the bandwith of the full matrix system.

The bandwidth calculated here corresponds to what is referred to as half-bandwidth. If the elements of the matrix are designated as a_ij, it corresponds to the maximum value of |i-j| for non-zero a_ij. As a result, the value also corresponds to the number of sub- or super-diagonals.

The bandwith can be calculated elementally as it corresponds to the maximal elemental bandwith (i.e. the maximal difference in global DOF index for every element).

We caluclate here the bandwith of the full global system.

Definition at line 499 of file AssemblyMapCG.cpp.

References m_fullSystemBandWidth, m_localToGlobalMap, Nektar::MultiRegions::AssemblyMap::m_numGlobalDirBndCoeffs, Nektar::MultiRegions::AssemblyMap::m_numLocalBndCoeffsPerPatch, Nektar::MultiRegions::AssemblyMap::m_numLocalCoeffs, Nektar::MultiRegions::AssemblyMap::m_numLocalIntCoeffsPerPatch, and Nektar::MultiRegions::AssemblyMap::m_numPatches.

Referenced by Nektar::MultiRegions::AssemblyMapCG1D::AssemblyMapCG1D(), Nektar::MultiRegions::AssemblyMapCG2D::AssemblyMapCG2D(), and Nektar::MultiRegions::AssemblyMapCG3D::AssemblyMapCG3D().

{
int i,j;
int cnt = 0;
int locSize;
int maxId;
int minId;
int bwidth = -1;
for(i = 0; i < m_numPatches; ++i)
{
maxId = -1;
minId = m_numLocalCoeffs+1;
for(j = 0; j < locSize; j++)
{
{
if(m_localToGlobalMap[cnt+j] > maxId)
{
maxId = m_localToGlobalMap[cnt+j];
}
if(m_localToGlobalMap[cnt+j] < minId)
{
minId = m_localToGlobalMap[cnt+j];
}
}
}
bwidth = (bwidth>(maxId-minId))?bwidth:(maxId-minId);
cnt+=locSize;
}
}
map<int, vector<ExtraDirDof> >& Nektar::MultiRegions::AssemblyMapCG::GetExtraDirDofs ( )
inline

Definition at line 82 of file AssemblyMapCG.h.

References m_extraDirDofs.

{
}
void Nektar::MultiRegions::AssemblyMapCG::SetUpUniversalC0ContMap ( const ExpList locExp,
const PeriodicMap perVerts = NullPeriodicMap,
const PeriodicMap perEdges = NullPeriodicMap,
const PeriodicMap perFaces = NullPeriodicMap 
)
protected

Sets up the global to universal mapping of degrees of freedom across processors.

Definition at line 152 of file AssemblyMapCG.cpp.

References Nektar::MultiRegions::DeterminePeriodicFaceOrient(), Nektar::MultiRegions::ExpList::GetCoeff_Offset(), Nektar::MultiRegions::ExpList::GetExp(), Nektar::StdRegions::StdExpansion::GetFaceOrient(), Gs::Init(), Nektar::MultiRegions::AssemblyMap::m_bndGsh, Nektar::MultiRegions::AssemblyMap::m_comm, Nektar::MultiRegions::AssemblyMap::m_globalToUniversalBndMap, Nektar::MultiRegions::AssemblyMap::m_globalToUniversalBndMapUnique, m_globalToUniversalMap, m_globalToUniversalMapUnique, Nektar::MultiRegions::AssemblyMap::m_gsh, m_localToGlobalMap, Nektar::MultiRegions::AssemblyMap::m_numGlobalBndCoeffs, Nektar::MultiRegions::AssemblyMap::m_numGlobalCoeffs, Nektar::LibUtilities::ReduceMax, Nektar::LibUtilities::ReduceSum, Gs::Unique(), and Vmath::Zero().

Referenced by Nektar::MultiRegions::AssemblyMapCG1D::SetUp1DExpansionC0ContMap(), Nektar::MultiRegions::AssemblyMapCG2D::SetUp2DExpansionC0ContMap(), and Nektar::MultiRegions::AssemblyMapCG3D::SetUp3DExpansionC0ContMap().

{
int nDim = 0;
int nVert = 0;
int nEdge = 0;
int nFace = 0;
int maxEdgeDof = 0;
int maxFaceDof = 0;
int maxIntDof = 0;
int dof = 0;
int cnt;
int i,j,k;
int meshVertId;
int meshEdgeId;
int meshFaceId;
int elementId;
int vGlobalId;
int maxBndGlobalId = 0;
Array<OneD, unsigned int> edgeInteriorMap;
Array<OneD, int> edgeInteriorSign;
Array<OneD, unsigned int> faceInteriorMap;
Array<OneD, int> faceInteriorSign;
Array<OneD, unsigned int> interiorMap;
PeriodicMap::const_iterator pIt;
const LocalRegions::ExpansionVector &locExpVector = *(locExp.GetExp());
LibUtilities::CommSharedPtr vCommRow = m_comm->GetRowComm();
m_globalToUniversalMap = Nektar::Array<OneD, int>(m_numGlobalCoeffs, -1);
m_globalToUniversalMapUnique = Nektar::Array<OneD, int>(m_numGlobalCoeffs, -1);
m_globalToUniversalBndMap = Nektar::Array<OneD, int>(m_numGlobalBndCoeffs, -1);
m_globalToUniversalBndMapUnique = Nektar::Array<OneD, int>(m_numGlobalBndCoeffs, -1);
// Loop over all the elements in the domain to gather mesh data
for(i = 0; i < locExpVector.size(); ++i)
{
locExpansion = locExpVector[i];
nVert += locExpansion->GetNverts();
nEdge += locExpansion->GetNedges();
nFace += locExpansion->GetNfaces();
// Loop over all edges (and vertices) of element i
for(j = 0; j < locExpansion->GetNedges(); ++j)
{
dof = locExpansion->GetEdgeNcoeffs(j)-2;
maxEdgeDof = (dof > maxEdgeDof ? dof : maxEdgeDof);
}
for(j = 0; j < locExpansion->GetNfaces(); ++j)
{
dof = locExpansion->GetFaceIntNcoeffs(j);
maxFaceDof = (dof > maxFaceDof ? dof : maxFaceDof);
}
locExpansion->GetInteriorMap(interiorMap);
dof = interiorMap.num_elements();
maxIntDof = (dof > maxIntDof ? dof : maxIntDof);
}
// Tell other processes about how many dof we have
vCommRow->AllReduce(nVert, LibUtilities::ReduceSum);
vCommRow->AllReduce(nEdge, LibUtilities::ReduceSum);
vCommRow->AllReduce(nFace, LibUtilities::ReduceSum);
vCommRow->AllReduce(maxEdgeDof, LibUtilities::ReduceMax);
vCommRow->AllReduce(maxFaceDof, LibUtilities::ReduceMax);
vCommRow->AllReduce(maxIntDof, LibUtilities::ReduceMax);
// Assemble global to universal mapping for this process
for(i = 0; i < locExpVector.size(); ++i)
{
locExpansion = locExpVector[i];
nDim = locExpansion->GetShapeDimension();
cnt = locExp.GetCoeff_Offset(i);
// Loop over all vertices of element i
for(j = 0; j < locExpansion->GetNverts(); ++j)
{
meshVertId = locExpansion->GetGeom()->GetVid(j);
vGlobalId = m_localToGlobalMap[cnt+locExpansion->GetVertexMap(j)];
pIt = perVerts.find(meshVertId);
if (pIt != perVerts.end())
{
for (k = 0; k < pIt->second.size(); ++k)
{
meshVertId = min(meshVertId, pIt->second[k].id);
}
}
m_globalToUniversalMap[vGlobalId] = meshVertId + 1;
m_globalToUniversalBndMap[vGlobalId]=m_globalToUniversalMap[vGlobalId];
maxBndGlobalId = (vGlobalId > maxBndGlobalId ? vGlobalId : maxBndGlobalId);
}
// Loop over all edges of element i
for(j = 0; j < locExpansion->GetNedges(); ++j)
{
meshEdgeId = locExpansion->GetGeom()->GetEid(j);
pIt = perEdges.find(meshEdgeId);
if (pIt != perEdges.end())
{
for (k = 0; k < pIt->second.size(); ++k)
{
meshEdgeId = min(meshEdgeId, pIt->second[k].id);
}
}
edgeOrient = locExpansion->GetGeom()->GetEorient(j);
locExpansion->GetEdgeInteriorMap(j,edgeOrient,edgeInteriorMap,edgeInteriorSign);
dof = locExpansion->GetEdgeNcoeffs(j)-2;
// Set the global DOF's for the interior modes of edge j
for(k = 0; k < dof; ++k)
{
vGlobalId = m_localToGlobalMap[cnt+edgeInteriorMap[k]];
= nVert + meshEdgeId * maxEdgeDof + k + 1;
m_globalToUniversalBndMap[vGlobalId]=m_globalToUniversalMap[vGlobalId];
maxBndGlobalId = (vGlobalId > maxBndGlobalId ? vGlobalId : maxBndGlobalId);
}
}
// Loop over all faces of element i
for(j = 0; j < locExpansion->GetNfaces(); ++j)
{
faceOrient = locExpansion->as<LocalRegions::Expansion3D>()
->GetGeom3D()->GetFaceOrient(j);
meshFaceId = locExpansion->GetGeom()->GetFid(j);
pIt = perFaces.find(meshFaceId);
if (pIt != perFaces.end())
{
if(meshFaceId == min(meshFaceId, pIt->second[0].id))
{
faceOrient = DeterminePeriodicFaceOrient(faceOrient,pIt->second[0].orient);
}
meshFaceId = min(meshFaceId, pIt->second[0].id);
}
locExpansion->GetFaceInteriorMap(j,faceOrient,faceInteriorMap,faceInteriorSign);
dof = locExpansion->GetFaceIntNcoeffs(j);
for(k = 0; k < dof; ++k)
{
vGlobalId = m_localToGlobalMap[cnt+faceInteriorMap[k]];
= nVert + nEdge*maxEdgeDof + meshFaceId * maxFaceDof
+ k + 1;
m_globalToUniversalBndMap[vGlobalId]=m_globalToUniversalMap[vGlobalId];
maxBndGlobalId = (vGlobalId > maxBndGlobalId ? vGlobalId : maxBndGlobalId);
}
}
// Add interior DOFs to complete universal numbering
locExpansion->GetInteriorMap(interiorMap);
dof = interiorMap.num_elements();
elementId = (locExpansion->GetGeom())->GetGlobalID();
for (k = 0; k < dof; ++k)
{
vGlobalId = m_localToGlobalMap[cnt+interiorMap[k]];
= nVert + nEdge*maxEdgeDof + nFace*maxFaceDof + elementId*maxIntDof + k + 1;
}
}
// Set up the GSLib universal assemble mapping
// Internal DOF do not participate in any data
// exchange, so we keep these set to the special GSLib id=0 so
// they are ignored.
Nektar::Array<OneD, long> tmp(m_numGlobalCoeffs);
Vmath::Zero(m_numGlobalCoeffs, tmp, 1);
Nektar::Array<OneD, long> tmp2(m_numGlobalBndCoeffs, tmp);
for (unsigned int i = 0; i < m_numGlobalBndCoeffs; ++i)
{
}
m_gsh = Gs::Init(tmp, vCommRow);
m_bndGsh = Gs::Init(tmp2, vCommRow);
Gs::Unique(tmp, vCommRow);
for (unsigned int i = 0; i < m_numGlobalCoeffs; ++i)
{
m_globalToUniversalMapUnique[i] = (tmp[i] >= 0 ? 1 : 0);
}
for (unsigned int i = 0; i < m_numGlobalBndCoeffs; ++i)
{
m_globalToUniversalBndMapUnique[i] = (tmp2[i] >= 0 ? 1 : 0);
}
}
void Nektar::MultiRegions::AssemblyMapCG::v_Assemble ( const Array< OneD, const NekDouble > &  loc,
Array< OneD, NekDouble > &  global 
) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 654 of file AssemblyMapCG.cpp.

References Vmath::Assmb(), m_localToGlobalMap, m_localToGlobalSign, Nektar::MultiRegions::AssemblyMap::m_numGlobalCoeffs, Nektar::MultiRegions::AssemblyMap::m_numLocalCoeffs, Nektar::MultiRegions::AssemblyMap::m_signChange, Nektar::MultiRegions::AssemblyMap::UniversalAssemble(), and Vmath::Zero().

{
Array<OneD, const NekDouble> local;
if(global.data() == loc.data())
{
local = Array<OneD, NekDouble>(local.num_elements(),local.data());
}
else
{
local = loc; // create reference
}
//ASSERTL1(loc.get() != global.get(),"Local and Global Arrays cannot be the same");
Vmath::Zero(m_numGlobalCoeffs, global.get(), 1);
{
Vmath::Assmb(m_numLocalCoeffs, m_localToGlobalSign.get(), local.get(), m_localToGlobalMap.get(), global.get());
}
else
{
Vmath::Assmb(m_numLocalCoeffs, local.get(), m_localToGlobalMap.get(), global.get());
}
}
void Nektar::MultiRegions::AssemblyMapCG::v_Assemble ( const NekVector< NekDouble > &  loc,
NekVector< NekDouble > &  global 
) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 682 of file AssemblyMapCG.cpp.

References Nektar::MultiRegions::AssemblyMap::Assemble(), and Nektar::NekVector< DataType >::GetPtr().

{
Assemble(loc.GetPtr(),global.GetPtr());
}
const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMapCG::v_GetExtraDirEdges ( )
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 751 of file AssemblyMapCG.cpp.

References m_extraDirEdges.

{
}
int Nektar::MultiRegions::AssemblyMapCG::v_GetFullSystemBandWidth ( ) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 711 of file AssemblyMapCG.cpp.

References m_fullSystemBandWidth.

{
}
int Nektar::MultiRegions::AssemblyMapCG::v_GetGlobalToUniversalMap ( const int  i) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 541 of file AssemblyMapCG.cpp.

References m_globalToUniversalMap.

{
}
const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMapCG::v_GetGlobalToUniversalMap ( void  )
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 558 of file AssemblyMapCG.cpp.

References m_globalToUniversalMap.

int Nektar::MultiRegions::AssemblyMapCG::v_GetGlobalToUniversalMapUnique ( const int  i) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 546 of file AssemblyMapCG.cpp.

References m_globalToUniversalMapUnique.

const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMapCG::v_GetGlobalToUniversalMapUnique ( void  )
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 564 of file AssemblyMapCG.cpp.

References m_globalToUniversalMapUnique.

int Nektar::MultiRegions::AssemblyMapCG::v_GetLocalToGlobalMap ( const int  i) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 536 of file AssemblyMapCG.cpp.

References m_localToGlobalMap.

{
return m_localToGlobalMap[i];
}
const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMapCG::v_GetLocalToGlobalMap ( void  )
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 552 of file AssemblyMapCG.cpp.

References m_localToGlobalMap.

{
}
NekDouble Nektar::MultiRegions::AssemblyMapCG::v_GetLocalToGlobalSign ( const int  i) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 569 of file AssemblyMapCG.cpp.

References m_localToGlobalSign, and Nektar::MultiRegions::AssemblyMap::m_signChange.

{
{
}
else
{
return 1.0;
}
}
const Array< OneD, NekDouble > & Nektar::MultiRegions::AssemblyMapCG::v_GetLocalToGlobalSign ( ) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 582 of file AssemblyMapCG.cpp.

References m_localToGlobalSign.

{
}
int Nektar::MultiRegions::AssemblyMapCG::v_GetNumDirEdges ( ) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 731 of file AssemblyMapCG.cpp.

References m_numDirEdges.

{
return m_numDirEdges;
}
int Nektar::MultiRegions::AssemblyMapCG::v_GetNumDirFaces ( ) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 736 of file AssemblyMapCG.cpp.

References m_numDirFaces.

{
return m_numDirFaces;
}
int Nektar::MultiRegions::AssemblyMapCG::v_GetNumNonDirEdgeModes ( ) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 721 of file AssemblyMapCG.cpp.

References m_numNonDirEdgeModes.

{
}
int Nektar::MultiRegions::AssemblyMapCG::v_GetNumNonDirEdges ( ) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 741 of file AssemblyMapCG.cpp.

References m_numNonDirEdges.

{
}
int Nektar::MultiRegions::AssemblyMapCG::v_GetNumNonDirFaceModes ( ) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 726 of file AssemblyMapCG.cpp.

References m_numNonDirFaceModes.

{
}
int Nektar::MultiRegions::AssemblyMapCG::v_GetNumNonDirFaces ( ) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 746 of file AssemblyMapCG.cpp.

References m_numNonDirFaces.

{
}
int Nektar::MultiRegions::AssemblyMapCG::v_GetNumNonDirVertexModes ( ) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 716 of file AssemblyMapCG.cpp.

References m_numNonDirVertexModes.

void Nektar::MultiRegions::AssemblyMapCG::v_GlobalToLocal ( const Array< OneD, const NekDouble > &  global,
Array< OneD, NekDouble > &  loc 
) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 622 of file AssemblyMapCG.cpp.

References Vmath::Gathr(), m_localToGlobalMap, m_localToGlobalSign, Nektar::MultiRegions::AssemblyMap::m_numLocalCoeffs, and Nektar::MultiRegions::AssemblyMap::m_signChange.

{
Array<OneD, const NekDouble> glo;
if(global.data() == loc.data())
{
glo = Array<OneD, NekDouble>(global.num_elements(),global.data());
}
else
{
glo = global; // create reference
}
{
}
else
{
Vmath::Gathr(m_numLocalCoeffs, glo.get(), m_localToGlobalMap.get(), loc.get());
}
}
void Nektar::MultiRegions::AssemblyMapCG::v_GlobalToLocal ( const NekVector< NekDouble > &  global,
NekVector< NekDouble > &  loc 
) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 647 of file AssemblyMapCG.cpp.

References Nektar::NekVector< DataType >::GetPtr(), and Nektar::MultiRegions::AssemblyMap::GlobalToLocal().

{
GlobalToLocal(global.GetPtr(),loc.GetPtr());
}
AssemblyMapSharedPtr Nektar::MultiRegions::AssemblyMapCG::v_LinearSpaceMap ( const ExpList locexp,
GlobalSysSolnType  solnType 
)
protectedvirtual

Construct an AssemblyMapCG object which corresponds to the linear space of the current object.

This function is used to create a linear-space assembly map, which is then used in the linear space preconditioner in the conjugate gradient solve.

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 358 of file AssemblyMapCG.cpp.

References Nektar::MultiRegions::eNull, Nektar::MultiRegions::ExpList::GetExp(), Gs::Init(), Nektar::MultiRegions::AssemblyMap::m_comm, m_globalToUniversalMap, m_localToGlobalMap, Nektar::MultiRegions::AssemblyMap::m_numGlobalCoeffs, Nektar::MultiRegions::AssemblyMap::m_numGlobalDirBndCoeffs, Nektar::MultiRegions::AssemblyMap::m_session, Gs::Unique(), and Vmath::Zero().

{
int i, j;
int nverts = 0;
const boost::shared_ptr<LocalRegions::ExpansionVector> exp
= locexp.GetExp();
int nelmts = exp->size();
// Get Default Map and turn off any searched values.
returnval->m_solnType = solnType;
returnval->m_preconType = eNull;
returnval->m_maxStaticCondLevel = 0;
returnval->m_signChange = false;
returnval->m_comm = m_comm;
// Count the number of vertices
for (i = 0; i < nelmts; ++i)
{
nverts += (*exp)[i]->GetNverts();
}
returnval->m_numLocalCoeffs = nverts;
returnval->m_localToGlobalMap = Array<OneD, int>(nverts, -1);
// Store original global ids in this map
returnval->m_localToGlobalBndMap = Array<OneD, int>(nverts, -1);
int cnt = 0;
int cnt1 = 0;
Array<OneD, int> GlobCoeffs(m_numGlobalCoeffs, -1);
// Set up local to global map;
for (i = 0; i < nelmts; ++i)
{
for (j = 0; j < (*exp)[i]->GetNverts(); ++j)
{
returnval->m_localToGlobalMap[cnt] =
returnval->m_localToGlobalBndMap[cnt] =
m_localToGlobalMap[cnt1 + (*exp)[i]->GetVertexMap(j,true)];
GlobCoeffs[returnval->m_localToGlobalMap[cnt]] = 1;
#if 1
// Set up numLocalDirBndCoeffs
if ((returnval->m_localToGlobalMap[cnt]) <
{
returnval->m_numLocalDirBndCoeffs++;
}
#endif
cnt++;
}
cnt1 += (*exp)[i]->GetNcoeffs();
}
cnt = 0;
// Reset global numbering and count number of dofs
for (i = 0; i < m_numGlobalCoeffs; ++i)
{
if (GlobCoeffs[i] != -1)
{
GlobCoeffs[i] = cnt++;
}
}
// Set up number of globalCoeffs;
returnval->m_numGlobalCoeffs = cnt;
// Set up number of global Dirichlet boundary coefficients
for (i = 0; i < m_numGlobalDirBndCoeffs; ++i)
{
if (GlobCoeffs[i] != -1)
{
returnval->m_numGlobalDirBndCoeffs++;
}
}
// Set up global to universal map
if (m_globalToUniversalMap.num_elements())
{
= m_session->GetComm()->GetRowComm();
int nglocoeffs = returnval->m_numGlobalCoeffs;
returnval->m_globalToUniversalMap
= Array<OneD, int> (nglocoeffs);
returnval->m_globalToUniversalMapUnique
= Array<OneD, int> (nglocoeffs);
// Reset local to global map and setup universal map
for (i = 0; i < nverts; ++i)
{
cnt = returnval->m_localToGlobalMap[i];
returnval->m_localToGlobalMap[i] = GlobCoeffs[cnt];
returnval->m_globalToUniversalMap[GlobCoeffs[cnt]] =
}
Nektar::Array<OneD, long> tmp(nglocoeffs);
Vmath::Zero(nglocoeffs, tmp, 1);
for (unsigned int i = 0; i < nglocoeffs; ++i)
{
tmp[i] = returnval->m_globalToUniversalMap[i];
}
returnval->m_gsh = Gs::Init(tmp, vCommRow);
Gs::Unique(tmp, vCommRow);
for (unsigned int i = 0; i < nglocoeffs; ++i)
{
returnval->m_globalToUniversalMapUnique[i]
= (tmp[i] >= 0 ? 1 : 0);
}
}
else // not sure this option is ever needed.
{
for (i = 0; i < nverts; ++i)
{
cnt = returnval->m_localToGlobalMap[i];
returnval->m_localToGlobalMap[i] = GlobCoeffs[cnt];
}
}
return returnval;
}
void Nektar::MultiRegions::AssemblyMapCG::v_LocalToGlobal ( const Array< OneD, const NekDouble > &  loc,
Array< OneD, NekDouble > &  global 
) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 587 of file AssemblyMapCG.cpp.

References Gs::Gather(), Gs::gs_max, Nektar::MultiRegions::AssemblyMap::m_gsh, m_localToGlobalMap, m_localToGlobalSign, Nektar::MultiRegions::AssemblyMap::m_numLocalCoeffs, Nektar::MultiRegions::AssemblyMap::m_signChange, and Vmath::Scatr().

{
Array<OneD, const NekDouble> local;
if(global.data() == loc.data())
{
local = Array<OneD, NekDouble>(loc.num_elements(),loc.data());
}
else
{
local = loc; // create reference
}
{
Vmath::Scatr(m_numLocalCoeffs, m_localToGlobalSign.get(), local.get(), m_localToGlobalMap.get(), global.get());
}
else
{
Vmath::Scatr(m_numLocalCoeffs, local.get(), m_localToGlobalMap.get(), global.get());
}
// ensure all values are unique by calling a max
}
void Nektar::MultiRegions::AssemblyMapCG::v_LocalToGlobal ( const NekVector< NekDouble > &  loc,
NekVector< NekDouble > &  global 
) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 615 of file AssemblyMapCG.cpp.

References Nektar::NekVector< DataType >::GetPtr(), and Nektar::MultiRegions::AssemblyMap::LocalToGlobal().

{
LocalToGlobal(loc.GetPtr(),global.GetPtr());
}
void Nektar::MultiRegions::AssemblyMapCG::v_UniversalAssemble ( Array< OneD, NekDouble > &  pGlobal) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 689 of file AssemblyMapCG.cpp.

References Gs::Gather(), Gs::gs_add, and Nektar::MultiRegions::AssemblyMap::m_gsh.

{
}
void Nektar::MultiRegions::AssemblyMapCG::v_UniversalAssemble ( NekVector< NekDouble > &  pGlobal) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 695 of file AssemblyMapCG.cpp.

References Nektar::NekVector< DataType >::GetPtr(), and Nektar::MultiRegions::AssemblyMap::UniversalAssemble().

{
}
void Nektar::MultiRegions::AssemblyMapCG::v_UniversalAssemble ( Array< OneD, NekDouble > &  pGlobal,
int  offset 
) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 701 of file AssemblyMapCG.cpp.

References Nektar::MultiRegions::AssemblyMap::UniversalAssemble(), and Vmath::Vcopy().

{
Array<OneD, NekDouble> tmp(offset);
Vmath::Vcopy(offset, pGlobal, 1, tmp, 1);
Vmath::Vcopy(offset, tmp, 1, pGlobal, 1);
}

Member Data Documentation

map<int, vector<ExtraDirDof> > Nektar::MultiRegions::AssemblyMapCG::m_extraDirDofs
protected

Map indicating degrees of freedom which are Dirichlet but whose value is stored on another processor.

Definition at line 119 of file AssemblyMapCG.h.

Referenced by GetExtraDirDofs(), Nektar::MultiRegions::AssemblyMapCG2D::SetUp2DExpansionC0ContMap(), and Nektar::MultiRegions::AssemblyMapCG3D::SetUp3DExpansionC0ContMap().

Array<OneD, int> Nektar::MultiRegions::AssemblyMapCG::m_extraDirEdges
protected

Extra dirichlet edges in parallel.

Definition at line 113 of file AssemblyMapCG.h.

Referenced by Nektar::MultiRegions::AssemblyMapCG3D::SetUp3DExpansionC0ContMap(), and v_GetExtraDirEdges().

int Nektar::MultiRegions::AssemblyMapCG::m_fullSystemBandWidth
protected

Bandwith of the full matrix system (no static condensation).

Definition at line 93 of file AssemblyMapCG.h.

Referenced by CalculateFullSystemBandWidth(), and v_GetFullSystemBandWidth().

Array<OneD,int> Nektar::MultiRegions::AssemblyMapCG::m_globalToUniversalMap
protected

Integer map of process coeffs to universal space.

Definition at line 95 of file AssemblyMapCG.h.

Referenced by SetUpUniversalC0ContMap(), v_GetGlobalToUniversalMap(), and v_LinearSpaceMap().

Array<OneD,int> Nektar::MultiRegions::AssemblyMapCG::m_globalToUniversalMapUnique
protected

Integer map of unique process coeffs to universal space (signed)

Definition at line 97 of file AssemblyMapCG.h.

Referenced by SetUpUniversalC0ContMap(), and v_GetGlobalToUniversalMapUnique().

Array<OneD,int> Nektar::MultiRegions::AssemblyMapCG::m_localToGlobalMap
protected

Integer map of local coeffs to global space.

Definition at line 89 of file AssemblyMapCG.h.

Referenced by CalculateFullSystemBandWidth(), Nektar::CoupledLocalToGlobalC0ContMap::CoupledLocalToGlobalC0ContMap(), Nektar::MultiRegions::AssemblyMapCG1D::SetUp1DExpansionC0ContMap(), Nektar::MultiRegions::AssemblyMapCG2D::SetUp2DExpansionC0ContMap(), Nektar::MultiRegions::AssemblyMapCG3D::SetUp3DExpansionC0ContMap(), SetUpUniversalC0ContMap(), v_Assemble(), v_GetLocalToGlobalMap(), v_GlobalToLocal(), v_LinearSpaceMap(), and v_LocalToGlobal().

Array<OneD,NekDouble> Nektar::MultiRegions::AssemblyMapCG::m_localToGlobalSign
protected

Integer sign of local coeffs to global space.

Definition at line 91 of file AssemblyMapCG.h.

Referenced by Nektar::CoupledLocalToGlobalC0ContMap::CoupledLocalToGlobalC0ContMap(), Nektar::MultiRegions::AssemblyMapCG2D::SetUp2DExpansionC0ContMap(), Nektar::MultiRegions::AssemblyMapCG3D::SetUp3DExpansionC0ContMap(), v_Assemble(), v_GetLocalToGlobalSign(), v_GlobalToLocal(), and v_LocalToGlobal().

int Nektar::MultiRegions::AssemblyMapCG::m_maxStaticCondLevel
protected

Maximum static condensation level.

Definition at line 116 of file AssemblyMapCG.h.

Referenced by AssemblyMapCG(), and Nektar::MultiRegions::AssemblyMapCG2D::SetUp2DExpansionC0ContMap().

int Nektar::MultiRegions::AssemblyMapCG::m_numDirEdges
protected

Number of Dirichlet edges.

Definition at line 105 of file AssemblyMapCG.h.

Referenced by Nektar::MultiRegions::AssemblyMapCG3D::SetUp3DExpansionC0ContMap(), and v_GetNumDirEdges().

int Nektar::MultiRegions::AssemblyMapCG::m_numDirFaces
protected

Number of Dirichlet faces.

Definition at line 107 of file AssemblyMapCG.h.

Referenced by Nektar::MultiRegions::AssemblyMapCG3D::SetUp3DExpansionC0ContMap(), and v_GetNumDirFaces().

int Nektar::MultiRegions::AssemblyMapCG::m_numNonDirEdgeModes
protected

Number of non Dirichlet edge modes.

Definition at line 101 of file AssemblyMapCG.h.

Referenced by Nektar::MultiRegions::AssemblyMapCG3D::SetUp3DExpansionC0ContMap(), and v_GetNumNonDirEdgeModes().

int Nektar::MultiRegions::AssemblyMapCG::m_numNonDirEdges
protected

Number of Dirichlet edges.

Definition at line 109 of file AssemblyMapCG.h.

Referenced by Nektar::MultiRegions::AssemblyMapCG2D::SetUp2DGraphC0ContMap(), Nektar::MultiRegions::AssemblyMapCG3D::SetUp3DExpansionC0ContMap(), and v_GetNumNonDirEdges().

int Nektar::MultiRegions::AssemblyMapCG::m_numNonDirFaceModes
protected

Number of non Dirichlet face modes.

Definition at line 103 of file AssemblyMapCG.h.

Referenced by Nektar::MultiRegions::AssemblyMapCG3D::SetUp3DExpansionC0ContMap(), and v_GetNumNonDirFaceModes().

int Nektar::MultiRegions::AssemblyMapCG::m_numNonDirFaces
protected

Number of Dirichlet faces.

Definition at line 111 of file AssemblyMapCG.h.

Referenced by Nektar::MultiRegions::AssemblyMapCG3D::SetUp3DExpansionC0ContMap(), and v_GetNumNonDirFaces().

int Nektar::MultiRegions::AssemblyMapCG::m_numNonDirVertexModes
protected

Number of non Dirichlet vertex modes.

Definition at line 99 of file AssemblyMapCG.h.

Referenced by Nektar::MultiRegions::AssemblyMapCG2D::SetUp2DGraphC0ContMap(), Nektar::MultiRegions::AssemblyMapCG3D::SetUp3DExpansionC0ContMap(), and v_GetNumNonDirVertexModes().