Nektar++
AssemblyMapDG.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File AssemblyMapDG.cpp
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: Local to Global Base Class mapping routines
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
37 #include <MultiRegions/ExpList.h>
38 #include <LocalRegions/SegExp.h>
39 #include <LocalRegions/TriExp.h>
40 #include <LocalRegions/QuadExp.h>
41 #include <LocalRegions/HexExp.h>
42 #include <LocalRegions/TetExp.h>
43 #include <LocalRegions/PrismExp.h>
44 #include <LocalRegions/PyrExp.h>
45 #include <LocalRegions/PointExp.h>
48 
49 #include <boost/core/ignore_unused.hpp>
50 #include <boost/config.hpp>
51 #include <boost/graph/adjacency_list.hpp>
52 #include <boost/graph/cuthill_mckee_ordering.hpp>
53 #include <boost/graph/properties.hpp>
54 #include <boost/graph/bandwidth.hpp>
55 
56 using namespace std;
57 
58 namespace Nektar
59 {
60  namespace MultiRegions
61  {
62  AssemblyMapDG::AssemblyMapDG():
63  m_numDirichletBndPhys(0)
64  {
65  }
66 
68  {
69  }
70 
74  const ExpListSharedPtr &trace,
75  const ExpList &locExp,
78  const PeriodicMap &periodicTrace,
79  const std::string variable):
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 
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;
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;
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 
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;
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  }
546 
547  /**
548  * Constructs a mapping between the process-local global numbering and
549  * a universal numbering of the trace space expansion. The universal
550  * numbering is defined by the mesh edge IDs to enforce consistency
551  * across processes.
552  *
553  * @param locExp List of local elemental expansions.
554  */
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]
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  }
732 
734  const ExpList &locExp,
735  const ExpListSharedPtr trace,
736  const PeriodicMap &perMap)
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  }
857 
859  Array<OneD, int> &toAlign,
861  int nquad1,
862  int nquad2)
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  }
937 
939  Array<OneD, NekDouble> &pGlobal) const
940  {
941  Gs::Gather(pGlobal, Gs::gs_add, m_traceGsh);
942  }
943 
945  {
946  return m_localToGlobalBndMap[i];
947  }
948 
950  {
951  return m_globalToUniversalBndMap[i];
952  }
953 
955  {
957  }
958 
960  {
961  return m_localToGlobalBndMap;
962  }
963 
965  {
967  }
968 
970  {
972  }
973 
975  const int i) const
976  {
977  return GetLocalToGlobalBndSign(i);
978  }
979 
982  Array<OneD, NekDouble>& global,
983  bool useComm ) const
984  {
985  boost::ignore_unused(useComm);
986  AssembleBnd(loc,global);
987  }
988 
990  const NekVector<NekDouble>& loc,
991  NekVector< NekDouble>& global,
992  bool useComm) const
993  {
994  boost::ignore_unused(useComm);
995  AssembleBnd(loc,global);
996  }
997 
999  const Array<OneD, const NekDouble>& global,
1000  Array<OneD, NekDouble>& loc) const
1001  {
1002  GlobalToLocalBnd(global,loc);
1003  }
1004 
1006  const NekVector<NekDouble>& global,
1007  NekVector< NekDouble>& loc) const
1008  {
1009  GlobalToLocalBnd(global,loc);
1010  }
1011 
1014  Array<OneD, NekDouble> &global) const
1015  {
1016  AssembleBnd(loc,global);
1017  }
1018 
1020  const NekVector<NekDouble>& loc,
1021  NekVector< NekDouble>& global) const
1022  {
1023  AssembleBnd(loc,global);
1024  }
1025 
1027  Array<OneD, NekDouble>& pGlobal) const
1028  {
1029  Gs::Gather(pGlobal, Gs::gs_add, m_gsh);
1030  }
1031 
1033  NekVector< NekDouble>& pGlobal) const
1034  {
1035  UniversalAssemble(pGlobal.GetPtr());
1036  }
1037 
1039  {
1040  return GetBndSystemBandWidth();
1041  }
1042 
1044  {
1045  return m_traceToUniversalMap[i];
1046  }
1047 
1049  {
1050  return m_traceToUniversalMapUnique[i];
1051  }
1052 
1054  {
1055  return m_numDirichletBndPhys;
1056  }
1057 
1060  {
1061  ASSERTL1(i >= 0 && i < m_elmtToTrace.num_elements(),
1062  "i is out of range");
1063  return m_elmtToTrace[i];
1064  }
1065 
1068  {
1069  return m_elmtToTrace;
1070  }
1071  } //namespace
1072 } // namespace
void GlobalToLocalBnd(const NekVector< NekDouble > &global, NekVector< NekDouble > &loc, int offset) const
AssemblyMapDG()
Default constructor.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
std::size_t hash_range(Iter first, Iter last)
Definition: HashUtils.hpp:69
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:163
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
int id
Geometry ID of entity.
void MultiLevelBisectionReordering(const BoostGraph &graph, Array< OneD, int > &perm, Array< OneD, int > &iperm, BottomUpSubStructuredGraphSharedPtr &substructgraph, std::set< int > partVerts, int mdswitch)
int GetBndSystemBandWidth() const
Returns the bandwidth of the boundary system.
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: AssemblyMap.h:307
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
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
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:167
Principle Modified Functions .
Definition: BasisType.h:48
Lagrange Polynomials using the Gauss points .
Definition: BasisType.h:55
STL namespace.
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)
void UniversalAssemble(Array< OneD, NekDouble > &pGlobal) const
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:2191
Array< OneD, int > m_traceToUniversalMap
Integer map of process trace space quadrature points to universal space.
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
Base class for constructing local to global mapping of degrees of freedom.
Definition: AssemblyMap.h:58
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
size_t m_hash
Hash for map.
Definition: AssemblyMap.h:310
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:103
Array< OneD, Array< OneD, LocalRegions::ExpansionSharedPtr > > & GetElmtToTrace()
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMapUnique()
void SetUpUniversalTraceMap(const ExpList &locExp, const ExpListSharedPtr trace, const PeriodicMap &perMap=NullPeriodicMap)
virtual void v_Assemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
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
virtual void v_GlobalToLocal(const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const
int m_numGlobalDirBndCoeffs
Number of Global Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:319
virtual void v_LocalToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, bool useComm) const
virtual int v_GetFullSystemBandWidth() const
std::shared_ptr< Geometry > GeometrySharedPtr
Definition: Geometry.h:53
virtual const Array< OneD, const int > & v_GetLocalToGlobalMap()
bool isLocal
Flag specifying if this entity is local to this partition.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
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, const NekDouble > GetLocalToGlobalBndSign() const
Retrieve the sign change for all local boundary modes.
void SetUpUniversalDGMap(const ExpList &locExp)
double NekDouble
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
The number of int dofs per patch.
Definition: AssemblyMap.h:391
std::map< int, std::vector< PeriodicEntity > > PeriodicMap
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
virtual ~AssemblyMapDG()
Destructor.
int m_numLocalDirBndCoeffs
Number of Local Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:317
Array< OneD, int > m_traceToUniversalMapUnique
Integer map of unique process trace space quadrature points to universal space (signed).
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
std::shared_ptr< PointExp > PointExpSharedPtr
Definition: PointExp.h:129
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
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
virtual const Array< OneD, NekDouble > & v_GetLocalToGlobalSign() const
void AssembleBnd(const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const
std::shared_ptr< SegExp > SegExpSharedPtr
Definition: SegExp.h:266
void UniversalTraceAssemble(Array< OneD, NekDouble > &pGlobal) const
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMap()
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
int m_numGlobalCoeffs
Total number of global coefficients.
Definition: AssemblyMap.h:343
Lagrange for SEM basis .
Definition: BasisType.h:54
#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
std::shared_ptr< BottomUpSubStructuredGraph > BottomUpSubStructuredGraphSharedPtr
std::shared_ptr< SessionReader > SessionReaderSharedPtr
void RealignTraceElement(Array< OneD, int > &toAlign, StdRegions::Orientation orient, int nquad1, int nquad2=0)
Array< OneD, DataType > & GetPtr()
Definition: NekVector.cpp:227
int GetNumDirichletBndPhys()
Return the number of boundary segments on which Dirichlet boundary conditions are imposed...
int m_numPatches
The number of patches (~elements) in the current level.
Definition: AssemblyMap.h:387
const SpatialDomains::PointGeomSharedPtr GetGeom() const
Definition: PointExp.h:110
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
virtual void v_UniversalAssemble(Array< OneD, NekDouble > &pGlobal) const