Nektar++
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
Nektar::MultiRegions::AssemblyMapDG Class Reference

#include <AssemblyMapDG.h>

Inheritance diagram for Nektar::MultiRegions::AssemblyMapDG:
[legend]

Public Member Functions

 AssemblyMapDG ()
 Default constructor. More...
 
 AssemblyMapDG (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &graph1D, const ExpListSharedPtr &trace, const ExpList &locExp, const Array< OneD, const MultiRegions::ExpListSharedPtr > &bndConstraint, const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > &bndCond, const PeriodicMap &periodicTrace, const std::string variable="DefaultVar")
 Constructor for trace map for one-dimensional expansion. More...
 
virtual ~AssemblyMapDG ()
 Destructor. More...
 
int GetNumDirichletBndPhys ()
 Return the number of boundary segments on which Dirichlet boundary conditions are imposed. More...
 
Array< OneD, LocalRegions::ExpansionSharedPtr > & GetElmtToTrace (const int i)
 
Array< OneD, Array< OneD, LocalRegions::ExpansionSharedPtr > > & GetElmtToTrace ()
 
int GetTraceToUniversalMap (int i)
 
int GetTraceToUniversalMapUnique (int i)
 
void UniversalTraceAssemble (Array< OneD, NekDouble > &pGlobal) const
 
- Public Member Functions inherited from Nektar::MultiRegions::AssemblyMap
 AssemblyMap ()
 Default constructor. More...
 
 AssemblyMap (const LibUtilities::SessionReaderSharedPtr &pSession, const std::string variable="DefaultVar")
 Constructor with a communicator. More...
 
 AssemblyMap (AssemblyMap *oldLevelMap, const BottomUpSubStructuredGraphSharedPtr &multiLevelGraph)
 Constructor for next level in multi-level static condensation. More...
 
virtual ~AssemblyMap ()
 Destructor. More...
 
LibUtilities::CommSharedPtr GetComm ()
 Retrieves the communicator. More...
 
size_t GetHash () const
 Retrieves the hash of this map. More...
 
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, bool useComm=true) const
 
void LocalToGlobal (const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, bool useComm=true) 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. More...
 
const Array< OneD, const int > & GetLocalToGlobalBndMap ()
 Retrieve the global indices of the local boundary modes. More...
 
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. More...
 
NekDouble GetLocalToGlobalBndSign (const int i) const
 Retrieve the sign change of a given local boundary mode. More...
 
Array< OneD, const NekDoubleGetLocalToGlobalBndSign () const
 Retrieve the sign change for all local boundary modes. More...
 
int GetBndCondCoeffsToGlobalCoeffsMap (const int i)
 Retrieves the global index corresponding to a boundary expansion mode. More...
 
const Array< OneD, const int > & GetBndCondCoeffsToGlobalCoeffsMap ()
 Retrieves the global indices corresponding to the boundary expansion modes. More...
 
NekDouble GetBndCondCoeffsToGlobalCoeffsSign (const int i)
 Returns the modal sign associated with a given boundary expansion mode. More...
 
int GetBndCondTraceToGlobalTraceMap (const int i)
 Returns the global index of the boundary trace giving the index on the boundary expansion. More...
 
const Array< OneD, const int > & GetBndCondTraceToGlobalTraceMap ()
 
int GetNumGlobalDirBndCoeffs () const
 Returns the number of global Dirichlet boundary coefficients. More...
 
int GetNumLocalDirBndCoeffs () const
 Returns the number of local Dirichlet boundary coefficients. More...
 
int GetNumGlobalBndCoeffs () const
 Returns the total number of global boundary coefficients. More...
 
int GetNumLocalBndCoeffs () const
 Returns the total number of local boundary coefficients. More...
 
int GetNumLocalCoeffs () const
 Returns the total number of local coefficients. More...
 
int GetNumGlobalCoeffs () const
 Returns the total number of global coefficients. More...
 
bool GetSingularSystem () const
 Retrieves if the system is singular (true) or not (false) More...
 
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, bool printHeader=true) const
 
const Array< OneD, const int > & GetExtraDirEdges ()
 
std::shared_ptr< AssemblyMapLinearSpaceMap (const ExpList &locexp, GlobalSysSolnType solnType)
 
int GetBndSystemBandWidth () const
 Returns the bandwidth of the boundary system. More...
 
int GetStaticCondLevel () const
 Returns the level of static condensation for this map. More...
 
int GetNumPatches () const
 Returns the number of patches in this static condensation level. More...
 
const Array< OneD, const unsigned int > & GetNumLocalBndCoeffsPerPatch ()
 Returns the number of local boundary coefficients in each patch. More...
 
const Array< OneD, const unsigned int > & GetNumLocalIntCoeffsPerPatch ()
 Returns the number of local interior coefficients in each patch. More...
 
const AssemblyMapSharedPtr GetNextLevelLocalToGlobalMap () const
 Returns the local to global mapping for the next level in the multi-level static condensation. More...
 
void SetNextLevelLocalToGlobalMap (AssemblyMapSharedPtr pNextLevelLocalToGlobalMap)
 
const PatchMapSharedPtrGetPatchMapFromPrevLevel (void) const
 Returns the patch map from the previous level of the multi-level static condensation. More...
 
bool AtLastLevel () const
 Returns true if this is the last level in the multi-level static condensation. More...
 
GlobalSysSolnType GetGlobalSysSolnType () const
 Returns the method of solving global systems. More...
 
PreconditionerType GetPreconType () const
 
NekDouble GetIterativeTolerance () const
 
int GetMaxIterations () const
 
int GetSuccessiveRHS () const
 
int GetLowestStaticCondLevel () const
 

Protected Member Functions

void SetUpUniversalDGMap (const ExpList &locExp)
 
void SetUpUniversalTraceMap (const ExpList &locExp, const ExpListSharedPtr trace, const PeriodicMap &perMap=NullPeriodicMap)
 
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 void v_LocalToGlobal (const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, bool useComm) const
 
virtual void v_LocalToGlobal (const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, bool useComm) 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 int v_GetFullSystemBandWidth () const
 
void RealignTraceElement (Array< OneD, int > &toAlign, StdRegions::Orientation orient, int nquad1, int nquad2=0)
 
- Protected Member Functions inherited from Nektar::MultiRegions::AssemblyMap
void CalculateBndSystemBandWidth ()
 Calculates the bandwidth of the boundary system. More...
 
void GlobalToLocalBndWithoutSign (const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc)
 

Protected Attributes

Gs::gs_datam_traceGsh
 
int m_numDirichletBndPhys
 Number of physical dirichlet boundary values in trace. More...
 
Array< OneD, Array< OneD, LocalRegions::ExpansionSharedPtr > > m_elmtToTrace
 list of edge expansions for a given element More...
 
Array< OneD, int > m_traceToUniversalMap
 Integer map of process trace space quadrature points to universal space. More...
 
Array< OneD, int > m_traceToUniversalMapUnique
 Integer map of unique process trace space quadrature points to universal space (signed). More...
 
- Protected Attributes inherited from Nektar::MultiRegions::AssemblyMap
LibUtilities::SessionReaderSharedPtr m_session
 Session object. More...
 
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
size_t m_hash
 Hash for map. More...
 
int m_numLocalBndCoeffs
 Number of local boundary coefficients. More...
 
int m_numGlobalBndCoeffs
 Total number of global boundary coefficients. More...
 
int m_numLocalDirBndCoeffs
 Number of Local Dirichlet Boundary Coefficients. More...
 
int m_numGlobalDirBndCoeffs
 Number of Global Dirichlet Boundary Coefficients. More...
 
bool m_systemSingular
 Flag indicating if the system is singular or not. More...
 
int m_numLocalCoeffs
 Total number of local coefficients. More...
 
int m_numGlobalCoeffs
 Total number of global coefficients. More...
 
bool m_signChange
 Flag indicating if modes require sign reversal. More...
 
Array< OneD, int > m_localToGlobalBndMap
 Integer map of local boundary coeffs to global space. More...
 
Array< OneD, NekDoublem_localToGlobalBndSign
 Integer sign of local boundary coeffs to global space. More...
 
Array< OneD, int > m_bndCondCoeffsToGlobalCoeffsMap
 Integer map of bnd cond coeffs to global coefficients. More...
 
Array< OneD, NekDoublem_bndCondCoeffsToGlobalCoeffsSign
 Integer map of bnd cond coeffs to global coefficients. More...
 
Array< OneD, int > m_bndCondTraceToGlobalTraceMap
 Integer map of bnd cond trace number to global trace number. More...
 
Array< OneD, int > m_globalToUniversalBndMap
 Integer map of process coeffs to universal space. More...
 
Array< OneD, int > m_globalToUniversalBndMapUnique
 Integer map of unique process coeffs to universal space (signed) More...
 
GlobalSysSolnType m_solnType
 The solution type of the global system. More...
 
int m_bndSystemBandWidth
 The bandwith of the global bnd system. More...
 
PreconditionerType m_preconType
 Type type of preconditioner to use in iterative solver. More...
 
int m_maxIterations
 Maximum iterations for iterative solver. More...
 
NekDouble m_iterativeTolerance
 Tolerance for iterative solver. More...
 
int m_successiveRHS
 sucessive RHS for iterative solver More...
 
Gs::gs_datam_gsh
 
Gs::gs_datam_bndGsh
 
int m_staticCondLevel
 The level of recursion in the case of multi-level static condensation. More...
 
int m_numPatches
 The number of patches (~elements) in the current level. More...
 
Array< OneD, unsigned int > m_numLocalBndCoeffsPerPatch
 The number of bnd dofs per patch. More...
 
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
 The number of int dofs per patch. More...
 
AssemblyMapSharedPtr m_nextLevelLocalToGlobalMap
 Map from the patches of the previous level to the patches of the current level. More...
 
int m_lowestStaticCondLevel
 Lowest static condensation level. More...
 

Detailed Description

Definition at line 52 of file AssemblyMapDG.h.

Constructor & Destructor Documentation

◆ AssemblyMapDG() [1/2]

Nektar::MultiRegions::AssemblyMapDG::AssemblyMapDG ( )

Default constructor.

Definition at line 62 of file AssemblyMapDG.cpp.

62  :
64  {
65  }
int m_numDirichletBndPhys
Number of physical dirichlet boundary values in trace.
Definition: AssemblyMapDG.h:96

◆ AssemblyMapDG() [2/2]

Nektar::MultiRegions::AssemblyMapDG::AssemblyMapDG ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr graph1D,
const ExpListSharedPtr trace,
const ExpList locExp,
const Array< OneD, const MultiRegions::ExpListSharedPtr > &  bndConstraint,
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > &  bndCond,
const PeriodicMap periodicTrace,
const std::string  variable = "DefaultVar" 
)

Constructor for trace map for one-dimensional expansion.

Definition at line 71 of file AssemblyMapDG.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, ASSERTL2, Nektar::MultiRegions::AssemblyMap::CalculateBndSystemBandWidth(), Nektar::MultiRegions::CuthillMckeeReordering(), Nektar::StdRegions::eDir1FwdDir1_Dir2FwdDir2, Nektar::MultiRegions::eDirectFullMatrix, Nektar::MultiRegions::eDirectMultiLevelStaticCond, Nektar::MultiRegions::eDirectStaticCond, Nektar::SpatialDomains::eDirichlet, Nektar::StdRegions::eForwards, Nektar::LibUtilities::eGauss_Lagrange, Nektar::LibUtilities::eGLL_Lagrange, Nektar::MultiRegions::eIterativeFull, Nektar::MultiRegions::eIterativeMultiLevelStaticCond, Nektar::MultiRegions::eIterativeStaticCond, Nektar::LibUtilities::eModified_A, Nektar::SpatialDomains::ePeriodic, Nektar::MultiRegions::ePETScFullMatrix, Nektar::MultiRegions::ePETScMultiLevelStaticCond, Nektar::MultiRegions::ePETScStaticCond, Nektar::MultiRegions::eXxtFullMatrix, Nektar::MultiRegions::eXxtMultiLevelStaticCond, Nektar::MultiRegions::eXxtStaticCond, Nektar::MultiRegions::ExpList::GetExp(), Nektar::hash_range(), Nektar::MultiRegions::AssemblyMap::m_bndCondCoeffsToGlobalCoeffsMap, Nektar::MultiRegions::AssemblyMap::m_bndCondTraceToGlobalTraceMap, m_elmtToTrace, Nektar::MultiRegions::AssemblyMap::m_hash, Nektar::MultiRegions::AssemblyMap::m_localToGlobalBndMap, Nektar::MultiRegions::AssemblyMap::m_localToGlobalBndSign, Nektar::MultiRegions::AssemblyMap::m_lowestStaticCondLevel, Nektar::MultiRegions::AssemblyMap::m_nextLevelLocalToGlobalMap, m_numDirichletBndPhys, Nektar::MultiRegions::AssemblyMap::m_numGlobalBndCoeffs, Nektar::MultiRegions::AssemblyMap::m_numGlobalCoeffs, Nektar::MultiRegions::AssemblyMap::m_numGlobalDirBndCoeffs, Nektar::MultiRegions::AssemblyMap::m_numLocalBndCoeffs, Nektar::MultiRegions::AssemblyMap::m_numLocalBndCoeffsPerPatch, Nektar::MultiRegions::AssemblyMap::m_numLocalCoeffs, Nektar::MultiRegions::AssemblyMap::m_numLocalDirBndCoeffs, Nektar::MultiRegions::AssemblyMap::m_numLocalIntCoeffsPerPatch, Nektar::MultiRegions::AssemblyMap::m_numPatches, Nektar::MultiRegions::AssemblyMap::m_signChange, Nektar::MultiRegions::AssemblyMap::m_solnType, Nektar::MultiRegions::AssemblyMap::m_staticCondLevel, Nektar::MultiRegions::MultiLevelBisectionReordering(), Nektar::MultiRegions::NoReordering(), SetUpUniversalDGMap(), and SetUpUniversalTraceMap().

79  :
80  AssemblyMap(pSession,variable)
81  {
82  boost::ignore_unused(graph);
83 
84  int i, j, k, cnt, eid, id, id1, gid;
85  int order_e = 0;
86  int nTraceExp = trace->GetExpSize();
87  int nbnd = bndCondExp.num_elements();
88 
92 
93  const LocalRegions::ExpansionVector expList = *(locExp.GetExp());
94  int nel = expList.size();
95 
96  map<int, int> meshTraceId;
97 
98  m_signChange = true;
99 
100  // determine mapping from geometry edges to trace
101  for(i = 0; i < nTraceExp; ++i)
102  {
103  meshTraceId[trace->GetExp(i)->GetGeom()->GetGlobalID()] = i;
104  }
105 
106  // Count total number of trace elements
107  cnt = 0;
108  for(i = 0; i < nel; ++i)
109  {
110  cnt += expList[i]->GetNtrace();
111  }
112 
113  Array<OneD, LocalRegions::ExpansionSharedPtr> tracemap(cnt);
114  m_elmtToTrace = Array<
115  OneD, Array<OneD, LocalRegions::ExpansionSharedPtr> >(nel);
116 
117  // set up trace expansions links;
118  for(cnt = i = 0; i < nel; ++i)
119  {
120  m_elmtToTrace[i] = tracemap + cnt;
121 
122  for(j = 0; j < expList[i]->GetNtrace(); ++j)
123  {
124  id = expList[i]->GetGeom()->GetTid(j);
125 
126  if(meshTraceId.count(id) > 0)
127  {
128  m_elmtToTrace[i][j] =
129  trace->GetExp(meshTraceId.find(id)->second);
130  }
131  else
132  {
133  ASSERTL0(false, "Failed to find trace map");
134  }
135  }
136 
137  cnt += expList[i]->GetNtrace();
138  }
139 
140  // Set up boundary mapping
141  cnt = 0;
142  for(i = 0; i < nbnd; ++i)
143  {
144  if (bndCond[i]->GetBoundaryConditionType() ==
146  {
147  continue;
148  }
149  cnt += bndCondExp[i]->GetExpSize();
150  }
151 
152  set<int> dirTrace;
153 
156 
157  for(i = 0; i < bndCondExp.num_elements(); ++i)
158  {
159  for(j = 0; j < bndCondExp[i]->GetExpSize(); ++j)
160  {
161  bndExp = bndCondExp[i]->GetExp(j);
162  traceGeom = bndExp->GetGeom();
163  id = traceGeom->GetGlobalID();
164 
165  if(bndCond[i]->GetBoundaryConditionType() ==
167  {
168  m_numLocalDirBndCoeffs += bndExp->GetNcoeffs();
169  m_numDirichletBndPhys += bndExp->GetTotPoints();
170  dirTrace.insert(id);
171  }
172  }
173  }
174 
175  // Set up integer mapping array and sign change for each degree of
176  // freedom + initialise some more data members.
177  m_staticCondLevel = 0;
179  m_numPatches = nel;
180  m_numLocalBndCoeffsPerPatch = Array<OneD, unsigned int>(nel);
181  m_numLocalIntCoeffsPerPatch = Array<OneD, unsigned int>(nel);
182 
183  int nbndry = 0;
184  for(i = 0; i < nel; ++i) // count number of elements in array
185  {
186  eid = i;
187  nbndry += expList[eid]->NumDGBndryCoeffs();
188  m_numLocalIntCoeffsPerPatch[i] = 0;
190  (unsigned int) expList[eid]->NumDGBndryCoeffs();
191  }
192 
194  m_numLocalBndCoeffs = nbndry;
195  m_numLocalCoeffs = nbndry;
196  m_localToGlobalBndMap = Array<OneD, int> (nbndry);
197  m_localToGlobalBndSign = Array<OneD, NekDouble> (nbndry,1);
198 
199  // Set up array for potential mesh optimsation
200  Array<OneD,int> traceElmtGid(nTraceExp, -1);
201  int nDir = 0;
202  cnt = 0;
203 
204  // We are now going to construct a graph of the mesh which can be
205  // reordered depending on the type of solver we would like to use.
206  typedef boost::adjacency_list<
207  boost::setS, boost::vecS, boost::undirectedS> BoostGraph;
208 
209  BoostGraph boostGraphObj;
210  int trace_id, trace_id1;
211  int dirOffset = 0;
212 
213  // make trace trace renumbering map where first solved trace starts
214  // at 0 so we can set up graph.
215  for(i = 0; i < nTraceExp; ++i)
216  {
217  id = trace->GetExp(i)->GetGeom()->GetGlobalID();
218 
219  if (dirTrace.count(id) == 0)
220  {
221  // Initial put in element ordering (starting from zero) into
222  // traceElmtGid
223  boost::add_vertex(boostGraphObj);
224  traceElmtGid[i] = cnt++;
225  }
226  else
227  {
228  // Use existing offset for Dirichlet edges
229  traceElmtGid[i] = dirOffset;
230  dirOffset += trace->GetExp(i)->GetNcoeffs();
231  nDir++;
232  }
233  }
234 
235  // Set up boost Graph
236  for(i = 0; i < nel; ++i)
237  {
238  eid = i;
239 
240  for(j = 0; j < expList[eid]->GetNtrace(); ++j)
241  {
242  // Add trace to boost graph for non-Dirichlet Boundary
243  traceGeom = m_elmtToTrace[eid][j]->GetGeom();
244  id = traceGeom->GetGlobalID();
245  trace_id = meshTraceId.find(id)->second;
246 
247  if(dirTrace.count(id) == 0)
248  {
249  for(k = j+1; k < expList[eid]->GetNtrace(); ++k)
250  {
251  traceGeom = m_elmtToTrace[eid][k]->GetGeom();
252  id1 = traceGeom->GetGlobalID();
253  trace_id1 = meshTraceId.find(id1)->second;
254 
255  if(dirTrace.count(id1) == 0)
256  {
257  boost::add_edge((size_t)traceElmtGid[trace_id],
258  (size_t)traceElmtGid[trace_id1],
259  boostGraphObj);
260  }
261  }
262  }
263  }
264  }
265 
266  int nGraphVerts = nTraceExp - nDir;
267  Array<OneD, int> perm (nGraphVerts);
268  Array<OneD, int> iperm(nGraphVerts);
269  Array<OneD, int> vwgts(nGraphVerts);
271 
272  for(i = 0; i < nGraphVerts; ++i)
273  {
274  vwgts[i] = trace->GetExp(i+nDir)->GetNcoeffs();
275  }
276 
277  if(nGraphVerts)
278  {
279  switch(m_solnType)
280  {
281  case eDirectFullMatrix:
282  case eIterativeFull:
284  case eXxtFullMatrix:
285  case eXxtStaticCond:
286  case ePETScFullMatrix:
287  case ePETScStaticCond:
288  {
289  NoReordering(boostGraphObj,perm,iperm);
290  break;
291  }
292  case eDirectStaticCond:
293  {
294  CuthillMckeeReordering(boostGraphObj,perm,iperm);
295  break;
296  }
301  {
302  MultiLevelBisectionReordering(boostGraphObj,perm,iperm,
303  bottomUpGraph);
304  break;
305  }
306  default:
307  {
308  ASSERTL0(false,"Unrecognised solution type");
309  }
310  }
311  }
312 
313  // Recast the permutation so that it can be used as a map from old
314  // trace ID to new trace ID
316  for(i = 0; i < nTraceExp - nDir; ++i)
317  {
318  traceElmtGid[perm[i]+nDir] = cnt;
319  cnt += trace->GetExp(perm[i]+nDir)->GetNcoeffs();
320  }
321 
322  // Now have trace edges Gid position
323 
324  cnt = 0;
325  for(i = 0; i < nel; ++i)
326  {
327  eid = i;
328  exp = expList[eid];
329 
330  for(j = 0; j < exp->GetNtrace(); ++j)
331  {
332  traceGeom = m_elmtToTrace[eid][j]->GetGeom();
333  id = traceGeom->GetGlobalID();
334  gid = traceElmtGid[meshTraceId.find(id)->second];
335 
336  const int nDim = expList[eid]->GetNumBases();
337 
338  if (nDim == 1)
339  {
340  order_e = 1;
341  m_localToGlobalBndMap[cnt] = gid;
342  }
343  else if (nDim == 2)
344  {
345  order_e = expList[eid]->GetEdgeNcoeffs(j);
346 
347  if(expList[eid]->GetEorient(j) == StdRegions::eForwards)
348  {
349  for(k = 0; k < order_e; ++k)
350  {
351  m_localToGlobalBndMap[k+cnt] = gid + k;
352  }
353  }
354  else
355  {
356  switch(m_elmtToTrace[eid][j]->GetBasisType(0))
357  {
359  {
360  // reverse vertex order
361  m_localToGlobalBndMap[cnt] = gid + 1;
362  m_localToGlobalBndMap[cnt+1] = gid;
363  for (k = 2; k < order_e; ++k)
364  {
365  m_localToGlobalBndMap[k+cnt] = gid + k;
366  }
367 
368  // negate odd modes
369  for(k = 3; k < order_e; k+=2)
370  {
371  m_localToGlobalBndSign[cnt+k] = -1.0;
372  }
373  break;
374  }
376  {
377  // reverse order
378  for(k = 0; k < order_e; ++k)
379  {
380  m_localToGlobalBndMap[cnt+order_e-k-1] = gid + k;
381  }
382  break;
383  }
385  {
386  // reverse order
387  for(k = 0; k < order_e; ++k)
388  {
389  m_localToGlobalBndMap[cnt+order_e-k-1] = gid + k;
390  }
391  break;
392  }
393  default:
394  {
395  ASSERTL0(false,"Boundary type not permitted");
396  }
397  }
398  }
399  }
400  else if (nDim == 3)
401  {
402  order_e = expList[eid]->GetFaceNcoeffs(j);
403 
404  std::map<int, int> orientMap;
405 
406  Array<OneD, unsigned int> elmMap1 (order_e);
407  Array<OneD, int> elmSign1(order_e);
408  Array<OneD, unsigned int> elmMap2 (order_e);
409  Array<OneD, int> elmSign2(order_e);
410 
411  StdRegions::Orientation fo = expList[eid]->GetForient(j);
412 
413  // Construct mapping which will permute global IDs
414  // according to face orientations.
415  expList[eid]->GetFaceToElementMap(j,fo,elmMap1,elmSign1);
416  expList[eid]->GetFaceToElementMap(
417  j,StdRegions::eDir1FwdDir1_Dir2FwdDir2,elmMap2,elmSign2);
418 
419  for (k = 0; k < elmMap1.num_elements(); ++k)
420  {
421  // Find the elemental co-efficient in the original
422  // mapping.
423  int idx = -1;
424  for (int l = 0; l < elmMap2.num_elements(); ++l)
425  {
426  if (elmMap1[k] == elmMap2[l])
427  {
428  idx = l;
429  break;
430  }
431  }
432 
433  ASSERTL2(idx != -1, "Problem with face to element map!");
434  orientMap[k] = idx;
435  }
436 
437  for(k = 0; k < order_e; ++k)
438  {
439  m_localToGlobalBndMap [k+cnt] = gid + orientMap[k];
440  m_localToGlobalBndSign[k+cnt] = elmSign2[orientMap[k]];
441  }
442  }
443 
444  cnt += order_e;
445  }
446  }
447 
448  // set up m_bndCondCoeffsToGlobalCoeffsMap to align with map
449  cnt = 0;
450  for(i = 0; i < nbnd; ++i)
451  {
452  if (bndCond[i]->GetBoundaryConditionType() ==
454  {
455  continue;
456  }
457  cnt += bndCondExp[i]->GetNcoeffs();
458  }
459 
460  m_bndCondCoeffsToGlobalCoeffsMap = Array<OneD, int>(cnt);
461 
462  // Number of boundary expansions
463  int nbndexp = 0, bndOffset, bndTotal = 0;
464  for(cnt = i = 0; i < nbnd; ++i)
465  {
466  if (bndCond[i]->GetBoundaryConditionType() ==
468  {
469  continue;
470  }
471 
472  for(j = 0; j < bndCondExp[i]->GetExpSize(); ++j)
473  {
474  bndExp = bndCondExp[i]->GetExp(j);
475  id = bndExp->GetGeom()->GetGlobalID();
476  gid = traceElmtGid[meshTraceId.find(id)->second];
477  bndOffset = bndCondExp[i]->GetCoeff_Offset(j) + bndTotal;
478 
479  // Since boundary information is defined to be aligned with
480  // the geometry just use forward/forward (both coordinate
481  // directions) defintiion for gid's
482  for(k = 0; k < bndExp->GetNcoeffs(); ++k)
483  {
484  m_bndCondCoeffsToGlobalCoeffsMap[bndOffset+k] = gid + k;
485  }
486  }
487 
488  nbndexp += bndCondExp[i]->GetExpSize();
489  bndTotal += bndCondExp[i]->GetNcoeffs();
490  }
491 
492  m_numGlobalBndCoeffs = trace->GetNcoeffs();
494 
496 
497  cnt = 0;
498  m_bndCondTraceToGlobalTraceMap = Array<OneD, int>(nbndexp);
499  for(i = 0; i < bndCondExp.num_elements(); ++i)
500  {
501  if (bndCond[i]->GetBoundaryConditionType() ==
503  {
504  continue;
505  }
506 
507  for(j = 0; j < bndCondExp[i]->GetExpSize(); ++j)
508  {
509  bndExp = bndCondExp[i]->GetExp(j);
510  traceGeom = bndExp->GetGeom();
511  id = traceGeom->GetGlobalID();
513  meshTraceId.find(id)->second;
514  }
515  }
516 
517  // Now set up mapping from global coefficients to universal.
518  ExpListSharedPtr tr = std::dynamic_pointer_cast<ExpList>(trace);
519  SetUpUniversalDGMap (locExp);
520  SetUpUniversalTraceMap(locExp, tr, periodicTrace);
521 
526  && nGraphVerts)
527  {
528  if (m_staticCondLevel < (bottomUpGraph->GetNlevels() - 1))
529  {
530  Array<OneD, int> vwgts_perm(nGraphVerts);
531 
532  for (int i = 0; i < nGraphVerts; i++)
533  {
534  vwgts_perm[i] = vwgts[perm[i]];
535  }
536 
537  bottomUpGraph->ExpandGraphWithVertexWeights(vwgts_perm);
539  AllocateSharedPtr(this, bottomUpGraph);
540  }
541  }
542 
544  m_localToGlobalBndMap.end());
545  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
std::size_t hash_range(Iter first, Iter last)
Definition: HashUtils.hpp:69
bool m_signChange
Flag indicating if modes require sign reversal.
Definition: AssemblyMap.h:346
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:315
void MultiLevelBisectionReordering(const BoostGraph &graph, Array< OneD, int > &perm, Array< OneD, int > &iperm, BottomUpSubStructuredGraphSharedPtr &substructgraph, std::set< int > partVerts, int mdswitch)
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
Array< OneD, int > m_bndCondTraceToGlobalTraceMap
Integer map of bnd cond trace number to global trace number.
Definition: AssemblyMap.h:357
Principle Modified Functions .
Definition: BasisType.h:48
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:55
int m_numLocalCoeffs
Total number of local coefficients.
Definition: AssemblyMap.h:332
void CuthillMckeeReordering(const BoostGraph &graph, Array< OneD, int > &perm, Array< OneD, int > &iperm)
std::vector< ExpansionSharedPtr > ExpansionVector
Definition: Expansion.h:67
AssemblyMapSharedPtr m_nextLevelLocalToGlobalMap
Map from the patches of the previous level to the patches of the current level.
Definition: AssemblyMap.h:396
size_t m_hash
Hash for map.
Definition: AssemblyMap.h:310
void SetUpUniversalTraceMap(const ExpList &locExp, const ExpListSharedPtr trace, const PeriodicMap &perMap=NullPeriodicMap)
Array< OneD, unsigned int > m_numLocalBndCoeffsPerPatch
The number of bnd dofs per patch.
Definition: AssemblyMap.h:389
GlobalSysSolnType m_solnType
The solution type of the global system.
Definition: AssemblyMap.h:364
int m_numGlobalDirBndCoeffs
Number of Global Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:319
std::shared_ptr< Geometry > GeometrySharedPtr
Definition: Geometry.h:53
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void SetUpUniversalDGMap(const ExpList &locExp)
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
The number of int dofs per patch.
Definition: AssemblyMap.h:391
int m_lowestStaticCondLevel
Lowest static condensation level.
Definition: AssemblyMap.h:398
void CalculateBndSystemBandWidth()
Calculates the bandwidth of the boundary system.
Array< OneD, int > m_localToGlobalBndMap
Integer map of local boundary coeffs to global space.
Definition: AssemblyMap.h:349
int m_numLocalDirBndCoeffs
Number of Local Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:317
Array< OneD, int > m_bndCondCoeffsToGlobalCoeffsMap
Integer map of bnd cond coeffs to global coefficients.
Definition: AssemblyMap.h:353
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:65
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:313
int m_staticCondLevel
The level of recursion in the case of multi-level static condensation.
Definition: AssemblyMap.h:385
Array< OneD, NekDouble > m_localToGlobalBndSign
Integer sign of local boundary coeffs to global space.
Definition: AssemblyMap.h:351
AssemblyMap()
Default constructor.
Definition: AssemblyMap.cpp:80
Array< OneD, Array< OneD, LocalRegions::ExpansionSharedPtr > > m_elmtToTrace
list of edge expansions for a given element
Definition: AssemblyMapDG.h:99
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:274
int m_numGlobalCoeffs
Total number of global coefficients.
Definition: AssemblyMap.h:343
Lagrange for SEM basis .
Definition: BasisType.h:54
std::shared_ptr< BottomUpSubStructuredGraph > BottomUpSubStructuredGraphSharedPtr
int m_numPatches
The number of patches (~elements) in the current level.
Definition: AssemblyMap.h:387
void NoReordering(const BoostGraph &graph, Array< OneD, int > &perm, Array< OneD, int > &iperm)
int m_numDirichletBndPhys
Number of physical dirichlet boundary values in trace.
Definition: AssemblyMapDG.h:96

◆ ~AssemblyMapDG()

Nektar::MultiRegions::AssemblyMapDG::~AssemblyMapDG ( )
virtual

Destructor.

Definition at line 67 of file AssemblyMapDG.cpp.

68  {
69  }

Member Function Documentation

◆ GetElmtToTrace() [1/2]

Array< OneD, LocalRegions::ExpansionSharedPtr > & Nektar::MultiRegions::AssemblyMapDG::GetElmtToTrace ( const int  i)

Definition at line 1059 of file AssemblyMapDG.cpp.

References ASSERTL1, and m_elmtToTrace.

1060  {
1061  ASSERTL1(i >= 0 && i < m_elmtToTrace.num_elements(),
1062  "i is out of range");
1063  return m_elmtToTrace[i];
1064  }
Array< OneD, Array< OneD, LocalRegions::ExpansionSharedPtr > > m_elmtToTrace
list of edge expansions for a given element
Definition: AssemblyMapDG.h:99
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250

◆ GetElmtToTrace() [2/2]

Array< OneD, Array< OneD, LocalRegions::ExpansionSharedPtr > > & Nektar::MultiRegions::AssemblyMapDG::GetElmtToTrace ( )

Definition at line 1067 of file AssemblyMapDG.cpp.

References m_elmtToTrace.

1068  {
1069  return m_elmtToTrace;
1070  }
Array< OneD, Array< OneD, LocalRegions::ExpansionSharedPtr > > m_elmtToTrace
list of edge expansions for a given element
Definition: AssemblyMapDG.h:99

◆ GetNumDirichletBndPhys()

int Nektar::MultiRegions::AssemblyMapDG::GetNumDirichletBndPhys ( )

Return the number of boundary segments on which Dirichlet boundary conditions are imposed.

Definition at line 1053 of file AssemblyMapDG.cpp.

References m_numDirichletBndPhys.

1054  {
1055  return m_numDirichletBndPhys;
1056  }
int m_numDirichletBndPhys
Number of physical dirichlet boundary values in trace.
Definition: AssemblyMapDG.h:96

◆ GetTraceToUniversalMap()

int Nektar::MultiRegions::AssemblyMapDG::GetTraceToUniversalMap ( int  i)

Definition at line 1043 of file AssemblyMapDG.cpp.

References m_traceToUniversalMap.

1044  {
1045  return m_traceToUniversalMap[i];
1046  }
Array< OneD, int > m_traceToUniversalMap
Integer map of process trace space quadrature points to universal space.

◆ GetTraceToUniversalMapUnique()

int Nektar::MultiRegions::AssemblyMapDG::GetTraceToUniversalMapUnique ( int  i)

Definition at line 1048 of file AssemblyMapDG.cpp.

References m_traceToUniversalMapUnique.

1049  {
1050  return m_traceToUniversalMapUnique[i];
1051  }
Array< OneD, int > m_traceToUniversalMapUnique
Integer map of unique process trace space quadrature points to universal space (signed).

◆ RealignTraceElement()

void Nektar::MultiRegions::AssemblyMapDG::RealignTraceElement ( Array< OneD, int > &  toAlign,
StdRegions::Orientation  orient,
int  nquad1,
int  nquad2 = 0 
)
protected

Definition at line 858 of file AssemblyMapDG.cpp.

References ASSERTL1, Nektar::StdRegions::eBackwards, Nektar::StdRegions::eDir1BwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1BwdDir1_Dir2FwdDir2, Nektar::StdRegions::eDir1BwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1BwdDir2_Dir2FwdDir1, Nektar::StdRegions::eDir1FwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1FwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1FwdDir2_Dir2FwdDir1, Nektar::StdRegions::eForwards, and Vmath::Vcopy().

Referenced by SetUpUniversalTraceMap().

863  {
864  if (orient == StdRegions::eBackwards)
865  {
866  ASSERTL1(nquad2 == 0, "nquad2 != 0 for reorienation");
867  for (int i = 0; i < nquad1/2; ++i)
868  {
869  swap(toAlign[i], toAlign[nquad1-1-i]);
870  }
871  }
872  else if (orient != StdRegions::eForwards)
873  {
874  ASSERTL1(nquad2 != 0, "nquad2 == 0 for reorienation");
875 
876  Array<OneD, int> tmp(nquad1*nquad2);
877 
878  // Copy transpose.
879  if (orient == StdRegions::eDir1FwdDir2_Dir2FwdDir1 ||
883  {
884  for (int i = 0; i < nquad2; ++i)
885  {
886  for (int j = 0; j < nquad1; ++j)
887  {
888  tmp[i*nquad1 + j] = toAlign[j*nquad2 + i];
889  }
890  }
891  }
892  else
893  {
894  for (int i = 0; i < nquad2; ++i)
895  {
896  for (int j = 0; j < nquad1; ++j)
897  {
898  tmp[i*nquad1 + j] = toAlign[i*nquad1 + j];
899  }
900  }
901  }
902 
903  if (orient == StdRegions::eDir1BwdDir1_Dir2FwdDir2 ||
907  {
908  // Reverse x direction
909  for (int i = 0; i < nquad2; ++i)
910  {
911  for (int j = 0; j < nquad1/2; ++j)
912  {
913  swap(tmp[i*nquad1 + j],
914  tmp[i*nquad1 + nquad1-j-1]);
915  }
916  }
917  }
918 
919  if (orient == StdRegions::eDir1FwdDir1_Dir2BwdDir2 ||
923  {
924  // Reverse y direction
925  for (int j = 0; j < nquad1; ++j)
926  {
927  for (int i = 0; i < nquad2/2; ++i)
928  {
929  swap(tmp[i*nquad1 + j],
930  tmp[(nquad2-i-1)*nquad1 + j]);
931  }
932  }
933  }
934  Vmath::Vcopy(nquad1*nquad2, tmp, 1, toAlign, 1);
935  }
936  }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064

◆ SetUpUniversalDGMap()

void Nektar::MultiRegions::AssemblyMapDG::SetUpUniversalDGMap ( const ExpList locExp)
protected

Constructs a mapping between the process-local global numbering and a universal numbering of the trace space expansion. The universal numbering is defined by the mesh edge IDs to enforce consistency across processes.

Parameters
locExpList of local elemental expansions.

Definition at line 555 of file AssemblyMapDG.cpp.

References ASSERTL2, Nektar::StdRegions::eDir1FwdDir1_Dir2FwdDir2, Nektar::StdRegions::eForwards, Nektar::MultiRegions::ExpList::GetExp(), Nektar::LocalRegions::Expansion::GetGeom(), Nektar::LocalRegions::PointExp::GetGeom(), Gs::Init(), Nektar::MultiRegions::AssemblyMap::m_bndGsh, Nektar::MultiRegions::AssemblyMap::m_comm, m_elmtToTrace, Nektar::MultiRegions::AssemblyMap::m_globalToUniversalBndMap, Nektar::MultiRegions::AssemblyMap::m_globalToUniversalBndMapUnique, Nektar::MultiRegions::AssemblyMap::m_gsh, Nektar::MultiRegions::AssemblyMap::m_localToGlobalBndMap, Nektar::MultiRegions::AssemblyMap::m_numGlobalBndCoeffs, Nektar::LibUtilities::ReduceMax, and Gs::Unique().

Referenced by AssemblyMapDG().

556  {
558  int eid = 0;
559  int cnt = 0;
560  int id = 0;
561  int order_e = 0;
562  int vGlobalId = 0;
563  int maxDof = 0;
564  int dof = 0;
565  int nDim = 0;
566  int i,j,k;
567 
568  const LocalRegions::ExpansionVector &locExpVector = *(locExp.GetExp());
569 
570  // Initialise the global to universal maps.
573 
574  // Loop over all the elements in the domain and compute max
575  // DOF. Reduce across all processes to get universal maximum.
576  for(i = 0; i < locExpVector.size(); ++i)
577  {
578  locExpansion = locExpVector[i];
579  nDim = locExpansion->GetShapeDimension();
580 
581  // Loop over all edges of element i
582  if (nDim == 1)
583  {
584  maxDof = (1 > maxDof ? 1 : maxDof);
585  }
586  else if (nDim == 2)
587  {
588  for (j = 0; j < locExpansion->GetNedges(); ++j)
589  {
590  dof = locExpansion->GetEdgeNcoeffs(j);
591  maxDof = (dof > maxDof ? dof : maxDof);
592  }
593  }
594  else if (nDim == 3)
595  {
596  for (j = 0; j < locExpansion->GetNfaces(); ++j)
597  {
598  dof = locExpansion->GetFaceNcoeffs(j);
599  maxDof = (dof > maxDof ? dof : maxDof);
600  }
601  }
602  }
603  m_comm->AllReduce(maxDof, LibUtilities::ReduceMax);
604 
605  // Now have trace edges Gid position
606  cnt = 0;
607  for(i = 0; i < locExpVector.size(); ++i)
608  {
609  eid = i;
610  locExpansion = locExpVector[eid];
611  nDim = locExpansion->GetShapeDimension();
612 
613  // Populate mapping for each edge of the element.
614  if (nDim == 1)
615  {
616  int nverts = locExpansion->GetNverts();
617  for(j = 0; j < nverts; ++j)
618  {
619  LocalRegions::PointExpSharedPtr locPointExp =
620  m_elmtToTrace[eid][j]->as<LocalRegions::PointExp>();
621  id = locPointExp->GetGeom()->GetGlobalID();
622  vGlobalId = m_localToGlobalBndMap[cnt+j];
623  m_globalToUniversalBndMap[vGlobalId]
624  = id * maxDof + j + 1;
625  }
626  cnt += nverts;
627  }
628  else if (nDim == 2)
629  {
630  for(j = 0; j < locExpansion->GetNedges(); ++j)
631  {
633  m_elmtToTrace[eid][j]->as<LocalRegions::SegExp>();
634 
635  id = locSegExp->GetGeom()->GetGlobalID();
636  order_e = locExpansion->GetEdgeNcoeffs(j);
637 
638  map<int,int> orientMap;
639  Array<OneD, unsigned int> map1(order_e), map2(order_e);
640  Array<OneD, int> sign1(order_e), sign2(order_e);
641 
642  locExpansion->GetEdgeToElementMap(j, StdRegions::eForwards, map1, sign1);
643  locExpansion->GetEdgeToElementMap(j, locExpansion->GetEorient(j), map2, sign2);
644 
645  for (k = 0; k < map1.num_elements(); ++k)
646  {
647  // Find the elemental co-efficient in the original
648  // mapping.
649  int idx = -1;
650  for (int l = 0; l < map2.num_elements(); ++l)
651  {
652  if (map1[k] == map2[l])
653  {
654  idx = l;
655  break;
656  }
657  }
658 
659  ASSERTL2(idx != -1, "Problem with face to element map!");
660  orientMap[k] = idx;
661  }
662 
663  for(k = 0; k < order_e; ++k)
664  {
665  vGlobalId = m_localToGlobalBndMap[k+cnt];
666  m_globalToUniversalBndMap[vGlobalId]
667  = id * maxDof + orientMap[k] + 1;
668  }
669  cnt += order_e;
670  }
671  }
672  else if (nDim == 3)
673  {
674  for(j = 0; j < locExpansion->GetNfaces(); ++j)
675  {
677  m_elmtToTrace[eid][j]
678  ->as<LocalRegions::Expansion2D>();
679 
680  id = locFaceExp->GetGeom()->GetGlobalID();
681  order_e = locExpansion->GetFaceNcoeffs(j);
682 
683  map<int,int> orientMap;
684  Array<OneD, unsigned int> map1(order_e), map2(order_e);
685  Array<OneD, int> sign1(order_e), sign2(order_e);
686 
687  locExpansion->GetFaceToElementMap(j, StdRegions::eDir1FwdDir1_Dir2FwdDir2, map1, sign1);
688  locExpansion->GetFaceToElementMap(j, locExpansion->GetForient(j), map2, sign2);
689 
690  for (k = 0; k < map1.num_elements(); ++k)
691  {
692  // Find the elemental co-efficient in the original
693  // mapping.
694  int idx = -1;
695  for (int l = 0; l < map2.num_elements(); ++l)
696  {
697  if (map1[k] == map2[l])
698  {
699  idx = l;
700  break;
701  }
702  }
703 
704  ASSERTL2(idx != -1, "Problem with face to element map!");
705  orientMap[k] = idx;
706  }
707 
708  for(k = 0; k < order_e; ++k)
709  {
710  vGlobalId = m_localToGlobalBndMap[k+cnt];
711  m_globalToUniversalBndMap[vGlobalId]
712  = id * maxDof + orientMap[k] + 1;
713  }
714  cnt += order_e;
715  }
716  }
717  }
718 
719  // Initialise GSlib and populate the unique map.
720  Array<OneD, long> tmp(m_globalToUniversalBndMap.num_elements());
721  for (i = 0; i < m_globalToUniversalBndMap.num_elements(); ++i)
722  {
723  tmp[i] = m_globalToUniversalBndMap[i];
724  }
725  m_bndGsh = m_gsh = Gs::Init(tmp, m_comm);
726  Gs::Unique(tmp, m_comm);
727  for (i = 0; i < m_globalToUniversalBndMap.num_elements(); ++i)
728  {
729  m_globalToUniversalBndMapUnique[i] = (tmp[i] >= 0 ? 1 : 0);
730  }
731  }
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:315
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: AssemblyMap.h:307
std::vector< ExpansionSharedPtr > ExpansionVector
Definition: Expansion.h:67
static gs_data * Init(const Nektar::Array< OneD, long > pId, const LibUtilities::CommSharedPtr &pComm, bool verbose=true)
Initialise Gather-Scatter map.
Definition: GsLib.hpp:167
static void Unique(const Nektar::Array< OneD, long > pId, const LibUtilities::CommSharedPtr &pComm)
Updates pId to negate all-but-one references to each universal ID.
Definition: GsLib.hpp:202
Array< OneD, int > m_localToGlobalBndMap
Integer map of local boundary coeffs to global space.
Definition: AssemblyMap.h:349
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:65
std::shared_ptr< PointExp > PointExpSharedPtr
Definition: PointExp.h:129
std::shared_ptr< Expansion2D > Expansion2DSharedPtr
Definition: Expansion1D.h:47
Array< OneD, Array< OneD, LocalRegions::ExpansionSharedPtr > > m_elmtToTrace
list of edge expansions for a given element
Definition: AssemblyMapDG.h:99
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:274
std::shared_ptr< SegExp > SegExpSharedPtr
Definition: SegExp.h:266
Array< OneD, int > m_globalToUniversalBndMap
Integer map of process coeffs to universal space.
Definition: AssemblyMap.h:359
Array< OneD, int > m_globalToUniversalBndMapUnique
Integer map of unique process coeffs to universal space (signed)
Definition: AssemblyMap.h:361

◆ SetUpUniversalTraceMap()

void Nektar::MultiRegions::AssemblyMapDG::SetUpUniversalTraceMap ( const ExpList locExp,
const ExpListSharedPtr  trace,
const PeriodicMap perMap = NullPeriodicMap 
)
protected

Definition at line 733 of file AssemblyMapDG.cpp.

References Nektar::MultiRegions::ExpList::GetExp(), Nektar::MultiRegions::PeriodicEntity::id, Gs::Init(), Nektar::MultiRegions::PeriodicEntity::isLocal, Nektar::MultiRegions::AssemblyMap::m_comm, m_traceGsh, m_traceToUniversalMap, m_traceToUniversalMapUnique, RealignTraceElement(), Nektar::LibUtilities::ReduceMax, and Gs::Unique().

Referenced by AssemblyMapDG().

737  {
738  Array<OneD, int> tmp;
740  int i;
741  int maxQuad = 0, quad = 0, nDim = 0, eid = 0, offset = 0;
742 
743  const LocalRegions::ExpansionVector &locExpVector = *(locExp.GetExp());
744 
745  int nTracePhys = trace->GetTotPoints();
746 
747  // Initialise the trace to universal maps.
749  Nektar::Array<OneD, int>(nTracePhys, -1);
751  Nektar::Array<OneD, int>(nTracePhys, -1);
752 
753  // Assume that each element of the expansion is of the same
754  // dimension.
755  nDim = locExpVector[0]->GetShapeDimension();
756 
757  if (nDim == 1)
758  {
759  maxQuad = (1 > maxQuad ? 1 : maxQuad);
760  }
761  else
762  {
763  for (i = 0; i < trace->GetExpSize(); ++i)
764  {
765  quad = trace->GetExp(i)->GetTotPoints();
766  if (quad > maxQuad)
767  {
768  maxQuad = quad;
769  }
770  }
771  }
772  m_comm->AllReduce(maxQuad, LibUtilities::ReduceMax);
773 
774  if (nDim == 1)
775  {
776  for (int i = 0; i < trace->GetExpSize(); ++i)
777  {
778  eid = trace->GetExp(i)->GetGeom()->GetGlobalID();
779  offset = trace->GetPhys_Offset(i);
780 
781  // Check to see if this vert is periodic. If it is, then we
782  // need use the unique eid of the two points
783  auto it = perMap.find(eid);
784  if (perMap.count(eid) > 0)
785  {
786  PeriodicEntity ent = it->second[0];
787  if (ent.isLocal == false) // Not sure if true in 1D
788  {
789  eid = min(eid, ent.id);
790  }
791  }
792 
793  m_traceToUniversalMap[offset] = eid*maxQuad+1;
794  }
795  }
796  else
797  {
798  for (int i = 0; i < trace->GetExpSize(); ++i)
799  {
800  eid = trace->GetExp(i)->GetGeom()->GetGlobalID();
801  offset = trace->GetPhys_Offset(i);
802  quad = trace->GetExp(i)->GetTotPoints();
803 
804  // Check to see if this edge is periodic. If it is, then we
805  // need to reverse the trace order of one edge only in the
806  // universal map so that the data are reversed w.r.t each
807  // other. We do this by using the minimum of the two IDs.
808  auto it = perMap.find(eid);
809  bool realign = false;
810  if (perMap.count(eid) > 0)
811  {
812  PeriodicEntity ent = it->second[0];
813  if (ent.isLocal == false)
814  {
815  realign = eid == min(eid, ent.id);
816  eid = min(eid, ent.id);
817  }
818  }
819 
820  for (int j = 0; j < quad; ++j)
821  {
822  m_traceToUniversalMap[j+offset] = eid*maxQuad+j+1;
823  }
824 
825  if (realign)
826  {
827  if (nDim == 2)
828  {
830  tmp = m_traceToUniversalMap+offset,
831  it->second[0].orient, quad);
832  }
833  else
834  {
836  tmp = m_traceToUniversalMap+offset,
837  it->second[0].orient,
838  trace->GetExp(i)->GetNumPoints(0),
839  trace->GetExp(i)->GetNumPoints(1));
840  }
841  }
842  }
843  }
844 
845  Array<OneD, long> tmp2(nTracePhys);
846  for (int i = 0; i < nTracePhys; ++i)
847  {
848  tmp2[i] = m_traceToUniversalMap[i];
849  }
850  m_traceGsh = Gs::Init(tmp2, m_comm);
851  Gs::Unique(tmp2, m_comm);
852  for (int i = 0; i < nTracePhys; ++i)
853  {
854  m_traceToUniversalMapUnique[i] = tmp2[i];
855  }
856  }
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: AssemblyMap.h:307
Array< OneD, int > m_traceToUniversalMap
Integer map of process trace space quadrature points to universal space.
std::vector< ExpansionSharedPtr > ExpansionVector
Definition: Expansion.h:67
static gs_data * Init(const Nektar::Array< OneD, long > pId, const LibUtilities::CommSharedPtr &pComm, bool verbose=true)
Initialise Gather-Scatter map.
Definition: GsLib.hpp:167
static void Unique(const Nektar::Array< OneD, long > pId, const LibUtilities::CommSharedPtr &pComm)
Updates pId to negate all-but-one references to each universal ID.
Definition: GsLib.hpp:202
Array< OneD, int > m_traceToUniversalMapUnique
Integer map of unique process trace space quadrature points to universal space (signed).
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:65
void RealignTraceElement(Array< OneD, int > &toAlign, StdRegions::Orientation orient, int nquad1, int nquad2=0)

◆ UniversalTraceAssemble()

void Nektar::MultiRegions::AssemblyMapDG::UniversalTraceAssemble ( Array< OneD, NekDouble > &  pGlobal) const

Definition at line 938 of file AssemblyMapDG.cpp.

References Gs::Gather(), Gs::gs_add, and m_traceGsh.

940  {
941  Gs::Gather(pGlobal, Gs::gs_add, m_traceGsh);
942  }
static void Gather(Nektar::Array< OneD, NekDouble > pU, gs_op pOp, gs_data *pGsh, Nektar::Array< OneD, NekDouble > pBuffer=NullNekDouble1DArray)
Performs a gather-scatter operation of the provided values.
Definition: GsLib.hpp:245

◆ v_Assemble() [1/2]

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

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 1012 of file AssemblyMapDG.cpp.

References Nektar::MultiRegions::AssemblyMap::AssembleBnd().

1015  {
1016  AssembleBnd(loc,global);
1017  }
void AssembleBnd(const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const

◆ v_Assemble() [2/2]

void Nektar::MultiRegions::AssemblyMapDG::v_Assemble ( const NekVector< NekDouble > &  loc,
NekVector< NekDouble > &  global 
) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 1019 of file AssemblyMapDG.cpp.

References Nektar::MultiRegions::AssemblyMap::AssembleBnd().

1022  {
1023  AssembleBnd(loc,global);
1024  }
void AssembleBnd(const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const

◆ v_GetFullSystemBandWidth()

int Nektar::MultiRegions::AssemblyMapDG::v_GetFullSystemBandWidth ( ) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 1038 of file AssemblyMapDG.cpp.

References Nektar::MultiRegions::AssemblyMap::GetBndSystemBandWidth().

1039  {
1040  return GetBndSystemBandWidth();
1041  }
int GetBndSystemBandWidth() const
Returns the bandwidth of the boundary system.

◆ v_GetGlobalToUniversalMap() [1/2]

int Nektar::MultiRegions::AssemblyMapDG::v_GetGlobalToUniversalMap ( const int  i) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 949 of file AssemblyMapDG.cpp.

References Nektar::MultiRegions::AssemblyMap::m_globalToUniversalBndMap.

950  {
951  return m_globalToUniversalBndMap[i];
952  }
Array< OneD, int > m_globalToUniversalBndMap
Integer map of process coeffs to universal space.
Definition: AssemblyMap.h:359

◆ v_GetGlobalToUniversalMap() [2/2]

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

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 964 of file AssemblyMapDG.cpp.

References Nektar::MultiRegions::AssemblyMap::m_globalToUniversalBndMap.

965  {
967  }
Array< OneD, int > m_globalToUniversalBndMap
Integer map of process coeffs to universal space.
Definition: AssemblyMap.h:359

◆ v_GetGlobalToUniversalMapUnique() [1/2]

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

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 954 of file AssemblyMapDG.cpp.

References Nektar::MultiRegions::AssemblyMap::m_globalToUniversalBndMapUnique.

955  {
957  }
Array< OneD, int > m_globalToUniversalBndMapUnique
Integer map of unique process coeffs to universal space (signed)
Definition: AssemblyMap.h:361

◆ v_GetGlobalToUniversalMapUnique() [2/2]

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

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 969 of file AssemblyMapDG.cpp.

References Nektar::MultiRegions::AssemblyMap::m_globalToUniversalBndMapUnique.

970  {
972  }
Array< OneD, int > m_globalToUniversalBndMapUnique
Integer map of unique process coeffs to universal space (signed)
Definition: AssemblyMap.h:361

◆ v_GetLocalToGlobalMap() [1/2]

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

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 944 of file AssemblyMapDG.cpp.

References Nektar::MultiRegions::AssemblyMap::m_localToGlobalBndMap.

945  {
946  return m_localToGlobalBndMap[i];
947  }
Array< OneD, int > m_localToGlobalBndMap
Integer map of local boundary coeffs to global space.
Definition: AssemblyMap.h:349

◆ v_GetLocalToGlobalMap() [2/2]

const Array< OneD, const int > & Nektar::MultiRegions::AssemblyMapDG::v_GetLocalToGlobalMap ( void  )
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 959 of file AssemblyMapDG.cpp.

References Nektar::MultiRegions::AssemblyMap::m_localToGlobalBndMap.

960  {
961  return m_localToGlobalBndMap;
962  }
Array< OneD, int > m_localToGlobalBndMap
Integer map of local boundary coeffs to global space.
Definition: AssemblyMap.h:349

◆ v_GetLocalToGlobalSign()

NekDouble Nektar::MultiRegions::AssemblyMapDG::v_GetLocalToGlobalSign ( const int  i) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 974 of file AssemblyMapDG.cpp.

References Nektar::MultiRegions::AssemblyMap::GetLocalToGlobalBndSign().

976  {
977  return GetLocalToGlobalBndSign(i);
978  }
Array< OneD, const NekDouble > GetLocalToGlobalBndSign() const
Retrieve the sign change for all local boundary modes.

◆ v_GlobalToLocal() [1/2]

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

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 998 of file AssemblyMapDG.cpp.

References Nektar::MultiRegions::AssemblyMap::GlobalToLocalBnd().

1001  {
1002  GlobalToLocalBnd(global,loc);
1003  }
void GlobalToLocalBnd(const NekVector< NekDouble > &global, NekVector< NekDouble > &loc, int offset) const

◆ v_GlobalToLocal() [2/2]

void Nektar::MultiRegions::AssemblyMapDG::v_GlobalToLocal ( const NekVector< NekDouble > &  global,
NekVector< NekDouble > &  loc 
) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 1005 of file AssemblyMapDG.cpp.

References Nektar::MultiRegions::AssemblyMap::GlobalToLocalBnd().

1008  {
1009  GlobalToLocalBnd(global,loc);
1010  }
void GlobalToLocalBnd(const NekVector< NekDouble > &global, NekVector< NekDouble > &loc, int offset) const

◆ v_LocalToGlobal() [1/2]

void Nektar::MultiRegions::AssemblyMapDG::v_LocalToGlobal ( const Array< OneD, const NekDouble > &  loc,
Array< OneD, NekDouble > &  global,
bool  useComm 
) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 980 of file AssemblyMapDG.cpp.

References Nektar::MultiRegions::AssemblyMap::AssembleBnd().

984  {
985  boost::ignore_unused(useComm);
986  AssembleBnd(loc,global);
987  }
void AssembleBnd(const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const

◆ v_LocalToGlobal() [2/2]

void Nektar::MultiRegions::AssemblyMapDG::v_LocalToGlobal ( const NekVector< NekDouble > &  loc,
NekVector< NekDouble > &  global,
bool  useComm 
) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 989 of file AssemblyMapDG.cpp.

References Nektar::MultiRegions::AssemblyMap::AssembleBnd().

993  {
994  boost::ignore_unused(useComm);
995  AssembleBnd(loc,global);
996  }
void AssembleBnd(const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const

◆ v_UniversalAssemble() [1/2]

void Nektar::MultiRegions::AssemblyMapDG::v_UniversalAssemble ( Array< OneD, NekDouble > &  pGlobal) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 1026 of file AssemblyMapDG.cpp.

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

1028  {
1029  Gs::Gather(pGlobal, Gs::gs_add, m_gsh);
1030  }
static void Gather(Nektar::Array< OneD, NekDouble > pU, gs_op pOp, gs_data *pGsh, Nektar::Array< OneD, NekDouble > pBuffer=NullNekDouble1DArray)
Performs a gather-scatter operation of the provided values.
Definition: GsLib.hpp:245

◆ v_UniversalAssemble() [2/2]

void Nektar::MultiRegions::AssemblyMapDG::v_UniversalAssemble ( NekVector< NekDouble > &  pGlobal) const
protectedvirtual

Reimplemented from Nektar::MultiRegions::AssemblyMap.

Definition at line 1032 of file AssemblyMapDG.cpp.

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

1034  {
1035  UniversalAssemble(pGlobal.GetPtr());
1036  }
void UniversalAssemble(Array< OneD, NekDouble > &pGlobal) const
Array< OneD, DataType > & GetPtr()
Definition: NekVector.cpp:227

Member Data Documentation

◆ m_elmtToTrace

Array<OneD, Array<OneD, LocalRegions::ExpansionSharedPtr> > Nektar::MultiRegions::AssemblyMapDG::m_elmtToTrace
protected

list of edge expansions for a given element

Definition at line 99 of file AssemblyMapDG.h.

Referenced by AssemblyMapDG(), GetElmtToTrace(), and SetUpUniversalDGMap().

◆ m_numDirichletBndPhys

int Nektar::MultiRegions::AssemblyMapDG::m_numDirichletBndPhys
protected

Number of physical dirichlet boundary values in trace.

Definition at line 96 of file AssemblyMapDG.h.

Referenced by AssemblyMapDG(), and GetNumDirichletBndPhys().

◆ m_traceGsh

Gs::gs_data* Nektar::MultiRegions::AssemblyMapDG::m_traceGsh
protected

Definition at line 93 of file AssemblyMapDG.h.

Referenced by SetUpUniversalTraceMap(), and UniversalTraceAssemble().

◆ m_traceToUniversalMap

Array<OneD,int> Nektar::MultiRegions::AssemblyMapDG::m_traceToUniversalMap
protected

Integer map of process trace space quadrature points to universal space.

Definition at line 102 of file AssemblyMapDG.h.

Referenced by GetTraceToUniversalMap(), and SetUpUniversalTraceMap().

◆ m_traceToUniversalMapUnique

Array<OneD,int> Nektar::MultiRegions::AssemblyMapDG::m_traceToUniversalMapUnique
protected

Integer map of unique process trace space quadrature points to universal space (signed).

Definition at line 105 of file AssemblyMapDG.h.

Referenced by GetTraceToUniversalMapUnique(), and SetUpUniversalTraceMap().