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 
38 #include <LocalRegions/PointExp.h>
39 #include <LocalRegions/SegExp.h>
42 #include <MultiRegions/ExpList.h>
43 
44 #include <boost/config.hpp>
45 #include <boost/core/ignore_unused.hpp>
46 #include <boost/graph/adjacency_list.hpp>
47 
48 using namespace std;
49 
50 namespace Nektar
51 {
52 namespace MultiRegions
53 {
54 AssemblyMapDG::AssemblyMapDG() : m_numDirichletBndPhys(0)
55 {
56 }
57 
59 {
60 }
61 
65  const ExpListSharedPtr &trace, const ExpList &locExp,
68  const PeriodicMap &periodicTrace, const std::string variable)
69  : AssemblyMap(pSession, locExp.GetComm(), variable)
70 {
71  boost::ignore_unused(graph);
72 
73  int i, j, k, cnt, id, id1, gid;
74  int order_e = 0;
75  int nTraceExp = trace->GetExpSize();
76  int nbnd = bndCondExp.size();
77 
81 
82  const LocalRegions::ExpansionVector expList = *(locExp.GetExp());
83  int nel = expList.size();
84 
85  map<int, int> meshTraceId;
86 
87  m_signChange = true;
88 
89  // determine mapping from geometry edges to trace
90  for (i = 0; i < nTraceExp; ++i)
91  {
92  meshTraceId[trace->GetExp(i)->GetGeom()->GetGlobalID()] = i;
93  }
94 
95  // Count total number of trace elements
96  cnt = 0;
97  for (i = 0; i < nel; ++i)
98  {
99  cnt += expList[i]->GetNtraces();
100  }
101 
103  m_elmtToTrace =
105 
106  // set up trace expansions links;
107  for (cnt = i = 0; i < nel; ++i)
108  {
109  m_elmtToTrace[i] = tracemap + cnt;
110  exp = expList[i];
111 
112  for (j = 0; j < exp->GetNtraces(); ++j)
113  {
114  id = exp->GetGeom()->GetTid(j);
115 
116  if (meshTraceId.count(id) > 0)
117  {
118  m_elmtToTrace[i][j] =
119  trace->GetExp(meshTraceId.find(id)->second);
120  }
121  else
122  {
123  ASSERTL0(false, "Failed to find trace map");
124  }
125  }
126 
127  cnt += exp->GetNtraces();
128  }
129 
130  // Set up boundary mapping
131  cnt = 0;
132  for (i = 0; i < nbnd; ++i)
133  {
134  if (bndCond[i]->GetBoundaryConditionType() == SpatialDomains::ePeriodic)
135  {
136  continue;
137  }
138  cnt += bndCondExp[i]->GetExpSize();
139  }
140 
141  set<int> dirTrace;
142 
145 
146  for (i = 0; i < bndCondExp.size(); ++i)
147  {
148  for (j = 0; j < bndCondExp[i]->GetExpSize(); ++j)
149  {
150  bndExp = bndCondExp[i]->GetExp(j);
151  traceGeom = bndExp->GetGeom();
152  id = traceGeom->GetGlobalID();
153 
154  if (bndCond[i]->GetBoundaryConditionType() ==
156  {
157  m_numLocalDirBndCoeffs += bndExp->GetNcoeffs();
158  m_numDirichletBndPhys += bndExp->GetTotPoints();
159  dirTrace.insert(id);
160  }
161  }
162  }
163 
164  // Set up integer mapping array and sign change for each degree of
165  // freedom + initialise some more data members.
166  m_staticCondLevel = 0;
168  m_numPatches = nel;
171 
172  int nbndry = 0;
173  for (i = 0; i < nel; ++i) // count number of elements in array
174  {
175  int BndCoeffs = expList[i]->NumDGBndryCoeffs();
176  nbndry += BndCoeffs;
178  m_numLocalBndCoeffsPerPatch[i] = (unsigned int)BndCoeffs;
179  }
180 
182  m_numLocalBndCoeffs = nbndry;
183  m_numLocalCoeffs = nbndry;
186 
187  // Set up array for potential mesh optimsation
188  Array<OneD, int> traceElmtGid(nTraceExp, -1);
189  int nDir = 0;
190  cnt = 0;
191 
192  // We are now going to construct a graph of the mesh which can be
193  // reordered depending on the type of solver we would like to use.
194  typedef boost::adjacency_list<boost::setS, boost::vecS, boost::undirectedS>
195  BoostGraph;
196 
197  BoostGraph boostGraphObj;
198  int trace_id, trace_id1;
199  int dirOffset = 0;
200 
201  // make trace renumbering map where first solved trace starts
202  // at 0 so we can set up graph.
203  for (i = 0; i < nTraceExp; ++i)
204  {
205  id = trace->GetExp(i)->GetGeom()->GetGlobalID();
206 
207  if (dirTrace.count(id) == 0)
208  {
209  // Initial put in element ordering (starting from zero) into
210  // traceElmtGid
211  boost::add_vertex(boostGraphObj);
212  traceElmtGid[i] = cnt++;
213  }
214  else
215  {
216  // Use existing offset for Dirichlet edges
217  traceElmtGid[i] = dirOffset;
218  dirOffset += trace->GetExp(i)->GetNcoeffs();
219  nDir++;
220  }
221  }
222 
223  // Set up boost Graph
224  for (i = 0; i < nel; ++i)
225  {
226  exp = expList[i];
227 
228  for (j = 0; j < exp->GetNtraces(); ++j)
229  {
230  // Add trace to boost graph for non-Dirichlet Boundary
231  traceGeom = m_elmtToTrace[i][j]->GetGeom();
232  id = traceGeom->GetGlobalID();
233  trace_id = meshTraceId.find(id)->second;
234 
235  if (dirTrace.count(id) == 0)
236  {
237  for (k = j + 1; k < exp->GetNtraces(); ++k)
238  {
239  traceGeom = m_elmtToTrace[i][k]->GetGeom();
240  id1 = traceGeom->GetGlobalID();
241  trace_id1 = meshTraceId.find(id1)->second;
242 
243  if (dirTrace.count(id1) == 0)
244  {
245  boost::add_edge((size_t)traceElmtGid[trace_id],
246  (size_t)traceElmtGid[trace_id1],
247  boostGraphObj);
248  }
249  }
250  }
251  }
252  }
253 
254  int nGraphVerts = nTraceExp - nDir;
255  Array<OneD, int> perm(nGraphVerts);
256  Array<OneD, int> iperm(nGraphVerts);
257  Array<OneD, int> vwgts(nGraphVerts);
259 
260  for (i = 0; i < nGraphVerts; ++i)
261  {
262  vwgts[i] = trace->GetExp(i + nDir)->GetNcoeffs();
263  }
264 
265  if (nGraphVerts)
266  {
267  switch (m_solnType)
268  {
269  case eDirectFullMatrix:
270  case eIterativeFull:
272  case eXxtFullMatrix:
273  case eXxtStaticCond:
274  case ePETScFullMatrix:
275  case ePETScStaticCond:
276  {
277  NoReordering(boostGraphObj, perm, iperm);
278  break;
279  }
280  case eDirectStaticCond:
281  {
282  CuthillMckeeReordering(boostGraphObj, perm, iperm);
283  break;
284  }
289  {
290  MultiLevelBisectionReordering(boostGraphObj, perm, iperm,
291  bottomUpGraph);
292  break;
293  }
294  default:
295  {
296  ASSERTL0(false, "Unrecognised solution type");
297  }
298  }
299  }
300 
301  // Recast the permutation so that it can be used as a map from old
302  // trace ID to new trace ID
304  for (i = 0; i < nTraceExp - nDir; ++i)
305  {
306  traceElmtGid[perm[i] + nDir] = cnt;
307  cnt += trace->GetExp(perm[i] + nDir)->GetNcoeffs();
308  }
309 
310  // Now have trace edges Gid position
311  cnt = 0;
312  for (i = 0; i < nel; ++i)
313  {
314  exp = expList[i];
315 
316  for (j = 0; j < exp->GetNtraces(); ++j)
317  {
318  traceGeom = m_elmtToTrace[i][j]->GetGeom();
319  id = traceGeom->GetGlobalID();
320  gid = traceElmtGid[meshTraceId.find(id)->second];
321 
322  const int nDim = exp->GetNumBases();
323 
324  if (nDim == 1)
325  {
326  order_e = 1;
327  m_localToGlobalBndMap[cnt] = gid;
328  }
329  else if (nDim == 2)
330  {
331  order_e = exp->GetTraceNcoeffs(j);
332 
333  if (exp->GetTraceOrient(j) == StdRegions::eForwards)
334  {
335  for (k = 0; k < order_e; ++k)
336  {
337  m_localToGlobalBndMap[k + cnt] = gid + k;
338  }
339  }
340  else
341  {
342  switch (m_elmtToTrace[i][j]->GetBasisType(0))
343  {
345  {
346  // reverse vertex order
347  m_localToGlobalBndMap[cnt] = gid + 1;
348  m_localToGlobalBndMap[cnt + 1] = gid;
349  for (k = 2; k < order_e; ++k)
350  {
351  m_localToGlobalBndMap[k + cnt] = gid + k;
352  }
353 
354  // negate odd modes
355  for (k = 3; k < order_e; k += 2)
356  {
357  m_localToGlobalBndSign[cnt + k] = -1.0;
358  }
359  break;
360  }
362  {
363  // reverse order
364  for (k = 0; k < order_e; ++k)
365  {
366  m_localToGlobalBndMap[cnt + order_e - k - 1] =
367  gid + k;
368  }
369  break;
370  }
372  {
373  // reverse order
374  for (k = 0; k < order_e; ++k)
375  {
376  m_localToGlobalBndMap[cnt + order_e - k - 1] =
377  gid + k;
378  }
379  break;
380  }
381  default:
382  {
383  ASSERTL0(false, "Boundary type not permitted");
384  }
385  }
386  }
387  }
388  else if (nDim == 3)
389  {
390  order_e = exp->GetTraceNcoeffs(j);
391 
392  std::map<int, int> orientMap;
393 
394  Array<OneD, unsigned int> elmMap1(order_e);
395  Array<OneD, int> elmSign1(order_e);
396  Array<OneD, unsigned int> elmMap2(order_e);
397  Array<OneD, int> elmSign2(order_e);
398 
399  StdRegions::Orientation fo = exp->GetTraceOrient(j);
400 
401  // Construct mapping which will permute global IDs
402  // according to face orientations.
403  exp->GetTraceToElementMap(j, elmMap1, elmSign1, fo);
404  exp->GetTraceToElementMap(j, elmMap2, elmSign2,
406 
407  for (k = 0; k < elmMap1.size(); ++k)
408  {
409  // Find the elemental co-efficient in the original
410  // mapping.
411  int idx = -1;
412  for (int l = 0; l < elmMap2.size(); ++l)
413  {
414  if (elmMap1[k] == elmMap2[l])
415  {
416  idx = l;
417  break;
418  }
419  }
420 
421  ASSERTL2(idx != -1, "Problem with face to element map!");
422  orientMap[k] = idx;
423  }
424 
425  for (k = 0; k < order_e; ++k)
426  {
427  m_localToGlobalBndMap[k + cnt] = gid + orientMap[k];
428  m_localToGlobalBndSign[k + cnt] = elmSign2[orientMap[k]];
429  }
430  }
431 
432  cnt += order_e;
433  }
434  }
435 
436  // set up identify map for lcoal to lcoal
438 
439  // local to bnd map is just a copy
440  for (i = 0; i < m_numLocalBndCoeffs; ++i)
441  {
442  m_localToLocalBndMap[i] = i;
443  }
444 
445  m_numGlobalBndCoeffs = trace->GetNcoeffs();
447 
449 
450  // set up m_bndCondCoeffsToLocalTraceMap
451  // Number of boundary expansions
452  int nbndexp = 0;
453  int bndTotal = 0;
454  int bndOffset;
455  int traceOffset;
456 
457  cnt = 0;
458  for (i = 0; i < nbnd; ++i)
459  {
460  if (bndCond[i]->GetBoundaryConditionType() == SpatialDomains::ePeriodic)
461  {
462  continue;
463  }
464  cnt += bndCondExp[i]->GetNcoeffs();
465  nbndexp += bndCondExp[i]->GetExpSize();
466  }
467 
470 
471  cnt = 0;
472  for (i = 0; i < bndCondExp.size(); ++i)
473  {
474  if (bndCond[i]->GetBoundaryConditionType() == SpatialDomains::ePeriodic)
475  {
476  continue;
477  }
478 
479  for (j = 0; j < bndCondExp[i]->GetExpSize(); ++j)
480  {
481  bndExp = bndCondExp[i]->GetExp(j);
482  id = bndExp->GetGeom()->GetGlobalID();
483 
484  int meshId = meshTraceId.find(id)->second;
485  m_bndCondIDToGlobalTraceID[cnt++] = meshId;
486 
487  // initialy set up map with global bnd location
488  // and then use the localToGlobalBndMap to invert
489  // since this is a one to one mapping on boundaries
490  traceOffset = traceElmtGid[meshId];
491  bndOffset = bndCondExp[i]->GetCoeff_Offset(j) + bndTotal;
492 
493  for (k = 0; k < bndExp->GetNcoeffs(); ++k)
494  {
495  m_bndCondCoeffsToLocalTraceMap[bndOffset + k] = traceOffset + k;
496  }
497  }
498  bndTotal += bndCondExp[i]->GetNcoeffs();
499  }
500 
501  // generate an inverse local to global bnd map;
502  map<int, int> invLocToGloMap;
503  for (i = 0; i < nbndry; ++i)
504  {
505  invLocToGloMap[m_localToGlobalBndMap[i]] = i;
506  }
507 
508  // reset bndCondCoeffToLocalTraceMap to hold local rather
509  // than global reference
510  for (i = 0; i < m_bndCondCoeffsToLocalTraceMap.size(); ++i)
511  {
513  invLocToGloMap[m_bndCondCoeffsToLocalTraceMap[i]];
514  }
515 
516  // Now set up mapping from global coefficients to universal.
517  ExpListSharedPtr tr = std::dynamic_pointer_cast<ExpList>(trace);
518  SetUpUniversalDGMap(locExp);
519 
522  locExp, tr, m_elmtToTrace, bndCondExp, bndCond, periodicTrace));
523 
528  nGraphVerts)
529  {
530  if (m_staticCondLevel < (bottomUpGraph->GetNlevels() - 1))
531  {
532  Array<OneD, int> vwgts_perm(nGraphVerts);
533 
534  for (int i = 0; i < nGraphVerts; i++)
535  {
536  vwgts_perm[i] = vwgts[perm[i]];
537  }
538 
539  bottomUpGraph->ExpandGraphWithVertexWeights(vwgts_perm);
542  bottomUpGraph);
543  }
544  }
545 
546  m_hash =
548 }
549 
550 /**
551  * Constructs a mapping between the process-local global numbering and
552  * a universal numbering of the trace space expansion. The universal
553  * numbering is defined by the mesh edge IDs to enforce consistency
554  * across processes.
555  *
556  * @param locExp List of local elemental expansions.
557  */
559 {
561  int cnt = 0;
562  int id = 0;
563  int order_e = 0;
564  int vGlobalId = 0;
565  int maxDof = 0;
566  int dof = 0;
567  int nDim = 0;
568  int i, j, k;
569 
570  const LocalRegions::ExpansionVector &locExpVector = *(locExp.GetExp());
571 
572  // Initialise the global to universal maps.
577 
578  // Loop over all the elements in the domain and compute max
579  // DOF. Reduce across all processes to get universal maximum.
580  for (i = 0; i < locExpVector.size(); ++i)
581  {
582  locExpansion = locExpVector[i];
583  nDim = locExpansion->GetShapeDimension();
584 
585  // Loop over all edges of element i
586  if (nDim == 1)
587  {
588  maxDof = (1 > maxDof ? 1 : maxDof);
589  }
590  else
591  {
592  for (j = 0; j < locExpansion->GetNtraces(); ++j)
593  {
594  dof = locExpansion->GetTraceNcoeffs(j);
595  maxDof = (dof > maxDof ? dof : maxDof);
596  }
597  }
598  }
599  m_comm->GetRowComm()->AllReduce(maxDof, LibUtilities::ReduceMax);
600 
601  // Now have trace edges Gid position
602  cnt = 0;
603  for (i = 0; i < locExpVector.size(); ++i)
604  {
605  locExpansion = locExpVector[i];
606  nDim = locExpansion->GetShapeDimension();
607 
608  // Populate mapping for each edge of the element.
609  if (nDim == 1)
610  {
611  int nverts = locExpansion->GetNverts();
612  for (j = 0; j < nverts; ++j)
613  {
614  LocalRegions::PointExpSharedPtr locPointExp =
616  id = locPointExp->GetGeom()->GetGlobalID();
617  vGlobalId = m_localToGlobalBndMap[cnt + j];
618  m_globalToUniversalBndMap[vGlobalId] = id * maxDof + j + 1;
619  }
620  cnt += nverts;
621  }
622  else if (nDim == 2)
623  {
624  for (j = 0; j < locExpansion->GetNtraces(); ++j)
625  {
627  m_elmtToTrace[i][j]->as<LocalRegions::SegExp>();
628 
629  id = locSegExp->GetGeom()->GetGlobalID();
630  order_e = locExpansion->GetTraceNcoeffs(j);
631 
632  map<int, int> orientMap;
633  Array<OneD, unsigned int> map1(order_e), map2(order_e);
634  Array<OneD, int> sign1(order_e), sign2(order_e);
635 
636  locExpansion->GetTraceToElementMap(j, map1, sign1,
638  locExpansion->GetTraceToElementMap(
639  j, map2, sign2, locExpansion->GetTraceOrient(j));
640 
641  for (k = 0; k < map1.size(); ++k)
642  {
643  // Find the elemental co-efficient in the original
644  // mapping.
645  int idx = -1;
646  for (int l = 0; l < map2.size(); ++l)
647  {
648  if (map1[k] == map2[l])
649  {
650  idx = l;
651  break;
652  }
653  }
654 
655  ASSERTL2(idx != -1, "Problem with face to"
656  " element map!");
657  orientMap[k] = idx;
658  }
659 
660  for (k = 0; k < order_e; ++k)
661  {
662  vGlobalId = m_localToGlobalBndMap[k + cnt];
663  m_globalToUniversalBndMap[vGlobalId] =
664  id * maxDof + orientMap[k] + 1;
665  }
666  cnt += order_e;
667  }
668  }
669  else if (nDim == 3) // This could likely be combined with nDim == 2
670  {
671  for (j = 0; j < locExpansion->GetNtraces(); ++j)
672  {
675 
676  id = locFaceExp->GetGeom()->GetGlobalID();
677  order_e = locExpansion->GetTraceNcoeffs(j);
678 
679  map<int, int> orientMap;
680  Array<OneD, unsigned int> map1(order_e), map2(order_e);
681  Array<OneD, int> sign1(order_e), sign2(order_e);
682 
683  locExpansion->GetTraceToElementMap(
684  j, map1, sign1, StdRegions::eDir1FwdDir1_Dir2FwdDir2);
685  locExpansion->GetTraceToElementMap(
686  j, map2, sign2, locExpansion->GetTraceOrient(j));
687 
688  for (k = 0; k < map1.size(); ++k)
689  {
690  // Find the elemental co-efficient in the original
691  // mapping.
692  int idx = -1;
693  for (int l = 0; l < map2.size(); ++l)
694  {
695  if (map1[k] == map2[l])
696  {
697  idx = l;
698  break;
699  }
700  }
701 
702  ASSERTL2(idx != -1, "Problem with face to "
703  "element map!");
704  orientMap[k] = idx;
705  }
706 
707  for (k = 0; k < order_e; ++k)
708  {
709  vGlobalId = m_localToGlobalBndMap[k + cnt];
710  m_globalToUniversalBndMap[vGlobalId] =
711  id * maxDof + orientMap[k] + 1;
712  }
713  cnt += order_e;
714  }
715  }
716  }
717 
718  // Initialise GSlib and populate the unique map.
720  for (i = 0; i < m_globalToUniversalBndMap.size(); ++i)
721  {
722  tmp[i] = m_globalToUniversalBndMap[i];
723  }
724  m_bndGsh = m_gsh = Gs::Init(tmp, m_comm->GetRowComm());
725  Gs::Unique(tmp, m_comm->GetRowComm());
726  for (i = 0; i < m_globalToUniversalBndMap.size(); ++i)
727  {
728  m_globalToUniversalBndMapUnique[i] = (tmp[i] >= 0 ? 1 : 0);
729  }
730 }
731 
734  int nquad1, int nquad2)
735 {
736  if (orient == StdRegions::eBackwards)
737  {
738  ASSERTL1(nquad2 == 0, "nquad2 != 0 for reorienation");
739  for (int i = 0; i < nquad1 / 2; ++i)
740  {
741  swap(toAlign[i], toAlign[nquad1 - 1 - i]);
742  }
743  }
744  else if (orient != StdRegions::eForwards)
745  {
746  ASSERTL1(nquad2 != 0, "nquad2 == 0 for reorienation");
747 
748  Array<OneD, int> tmp(nquad1 * nquad2);
749 
750  // Copy transpose.
751  if (orient == StdRegions::eDir1FwdDir2_Dir2FwdDir1 ||
755  {
756  for (int i = 0; i < nquad2; ++i)
757  {
758  for (int j = 0; j < nquad1; ++j)
759  {
760  tmp[i * nquad1 + j] = toAlign[j * nquad2 + i];
761  }
762  }
763  }
764  else
765  {
766  for (int i = 0; i < nquad2; ++i)
767  {
768  for (int j = 0; j < nquad1; ++j)
769  {
770  tmp[i * nquad1 + j] = toAlign[i * nquad1 + j];
771  }
772  }
773  }
774 
775  if (orient == StdRegions::eDir1BwdDir1_Dir2FwdDir2 ||
779  {
780  // Reverse x direction
781  for (int i = 0; i < nquad2; ++i)
782  {
783  for (int j = 0; j < nquad1 / 2; ++j)
784  {
785  swap(tmp[i * nquad1 + j], tmp[i * nquad1 + nquad1 - j - 1]);
786  }
787  }
788  }
789 
790  if (orient == StdRegions::eDir1FwdDir1_Dir2BwdDir2 ||
794  {
795  // Reverse y direction
796  for (int j = 0; j < nquad1; ++j)
797  {
798  for (int i = 0; i < nquad2 / 2; ++i)
799  {
800  swap(tmp[i * nquad1 + j],
801  tmp[(nquad2 - i - 1) * nquad1 + j]);
802  }
803  }
804  }
805  Vmath::Vcopy(nquad1 * nquad2, tmp, 1, toAlign, 1);
806  }
807 }
808 
810 {
811  return m_localToGlobalBndMap[i];
812 }
813 
815 {
816  return m_globalToUniversalBndMap[i];
817 }
818 
820 {
822 }
823 
825 {
826  return m_localToGlobalBndMap;
827 }
828 
830 {
832 }
833 
835 {
837 }
838 
840 {
841  return GetLocalToGlobalBndSign(i);
842 }
843 
845  Array<OneD, NekDouble> &global,
846  bool useComm) const
847 {
848  LocalBndToGlobal(loc, global, useComm);
849 }
850 
853 {
854  GlobalToLocalBnd(global, loc);
855 }
856 
858  NekVector<NekDouble> &loc) const
859 {
860  GlobalToLocalBnd(global, loc);
861 }
862 
864  Array<OneD, NekDouble> &global) const
865 {
866  AssembleBnd(loc, global);
867 }
868 
870  NekVector<NekDouble> &global) const
871 {
872  AssembleBnd(loc, global);
873 }
874 
876 {
877  Gs::Gather(pGlobal, Gs::gs_add, m_gsh);
878 }
879 
881 {
882  UniversalAssemble(pGlobal.GetPtr());
883 }
884 
886 {
887  return GetBndSystemBandWidth();
888 }
889 
891 {
892  return m_numDirichletBndPhys;
893 }
894 
896  const int i)
897 {
898  ASSERTL1(i >= 0 && i < m_elmtToTrace.size(), "i is out of range");
899  return m_elmtToTrace[i];
900 }
901 
904 {
905  return m_elmtToTrace;
906 }
907 
909 {
910  return m_assemblyComm;
911 }
912 
913 } // namespace MultiRegions
914 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:272
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:171
const SpatialDomains::PointGeomSharedPtr GetGeom() const
Definition: PointExp.h:113
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
Array< OneD, Array< OneD, LocalRegions::ExpansionSharedPtr > > m_elmtToTrace
list of edge expansions for a given element
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMap() override
virtual void v_LocalToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, bool useComm=false) const override
virtual void v_GlobalToLocal(const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const override
static void RealignTraceElement(Array< OneD, int > &toAlign, StdRegions::Orientation orient, int nquad1, int nquad2=0)
AssemblyCommDGSharedPtr m_assemblyComm
Array< OneD, Array< OneD, LocalRegions::ExpansionSharedPtr > > & GetElmtToTrace()
virtual const Array< OneD, const int > & v_GetLocalToGlobalMap() override
AssemblyCommDGSharedPtr GetAssemblyCommDG()
int GetNumDirichletBndPhys()
Return the number of boundary segments on which Dirichlet boundary conditions are imposed.
virtual void v_Assemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const override
void SetUpUniversalDGMap(const ExpList &locExp)
int m_numDirichletBndPhys
Number of physical dirichlet boundary values in trace.
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMapUnique() override
virtual void v_UniversalAssemble(Array< OneD, NekDouble > &pGlobal) const override
virtual int v_GetFullSystemBandWidth() const override
AssemblyMapDG()
Default constructor.
Base class for constructing local to global mapping of degrees of freedom.
Definition: AssemblyMap.h:58
int m_lowestStaticCondLevel
Lowest static condensation level.
Definition: AssemblyMap.h:437
GlobalSysSolnType m_solnType
The solution type of the global system.
Definition: AssemblyMap.h:398
int m_numLocalCoeffs
Total number of local coefficients.
Definition: AssemblyMap.h:360
Array< OneD, int > m_bndCondCoeffsToLocalTraceMap
Integer map of bnd cond coeff to local trace coeff.
Definition: AssemblyMap.h:389
bool m_signChange
Flag indicating if modes require sign reversal.
Definition: AssemblyMap.h:374
void AssembleBnd(const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const
int m_numGlobalCoeffs
Total number of global coefficients.
Definition: AssemblyMap.h:371
Array< OneD, int > m_globalToUniversalBndMap
Integer map of process coeffs to universal space.
Definition: AssemblyMap.h:393
void CalculateBndSystemBandWidth()
Calculates the bandwidth of the boundary system.
Array< OneD, NekDouble > m_localToGlobalBndSign
Integer sign of local boundary coeffs to global space.
Definition: AssemblyMap.h:379
Array< OneD, unsigned int > m_numLocalBndCoeffsPerPatch
The number of bnd dofs per patch.
Definition: AssemblyMap.h:428
virtual const Array< OneD, NekDouble > & v_GetLocalToGlobalSign() const
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:341
AssemblyMapSharedPtr m_nextLevelLocalToGlobalMap
Map from the patches of the previous level to the patches of the current level.
Definition: AssemblyMap.h:435
int m_staticCondLevel
The level of recursion in the case of multi-level static condensation.
Definition: AssemblyMap.h:424
int m_numLocalDirBndCoeffs
Number of Local Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:345
int m_numGlobalDirBndCoeffs
Number of Global Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:347
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
The number of int dofs per patch.
Definition: AssemblyMap.h:430
int GetBndSystemBandWidth() const
Returns the bandwidth of the boundary system.
Array< OneD, int > m_localToGlobalBndMap
Integer map of local coeffs to global Boundary Dofs.
Definition: AssemblyMap.h:377
Array< OneD, const NekDouble > GetLocalToGlobalBndSign() const
Retrieve the sign change for all local boundary modes.
Array< OneD, int > m_globalToUniversalBndMapUnique
Integer map of unique process coeffs to universal space (signed)
Definition: AssemblyMap.h:395
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: AssemblyMap.h:335
Array< OneD, int > m_localToLocalBndMap
Integer map of local boundary coeffs to local boundary system numbering.
Definition: AssemblyMap.h:381
void LocalBndToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, int offset, bool UseComm=true) const
void UniversalAssemble(Array< OneD, NekDouble > &pGlobal) const
int m_numPatches
The number of patches (~elements) in the current level.
Definition: AssemblyMap.h:426
void GlobalToLocalBnd(const NekVector< NekDouble > &global, NekVector< NekDouble > &loc, int offset) const
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:343
Array< OneD, int > m_bndCondIDToGlobalTraceID
Integer map of bnd cond trace number to global trace number.
Definition: AssemblyMap.h:391
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:103
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:2029
Array< OneD, DataType > & GetPtr()
Definition: NekVector.cpp:217
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:270
@ gs_add
Definition: GsLib.hpp:62
static gs_data * Init(const Nektar::Array< OneD, long > pId, const LibUtilities::CommSharedPtr &pComm, bool verbose=true)
Initialise Gather-Scatter map.
Definition: GsLib.hpp:192
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:227
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ eGauss_Lagrange
Lagrange Polynomials using the Gauss points.
Definition: BasisType.h:59
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:58
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50
std::shared_ptr< SegExp > SegExpSharedPtr
Definition: SegExp.h:251
std::shared_ptr< PointExp > PointExpSharedPtr
Definition: PointExp.h:135
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
std::shared_ptr< Expansion2D > Expansion2DSharedPtr
Definition: Expansion1D.h:47
std::vector< ExpansionSharedPtr > ExpansionVector
Definition: Expansion.h:70
void NoReordering(const BoostGraph &graph, Array< OneD, int > &perm, Array< OneD, int > &iperm)
std::shared_ptr< AssemblyCommDG > AssemblyCommDGSharedPtr
void CuthillMckeeReordering(const BoostGraph &graph, Array< OneD, int > &perm, Array< OneD, int > &iperm)
void MultiLevelBisectionReordering(const BoostGraph &graph, Array< OneD, int > &perm, Array< OneD, int > &iperm, BottomUpSubStructuredGraphSharedPtr &substructgraph, std::set< int > partVerts, int mdswitch)
std::shared_ptr< BottomUpSubStructuredGraph > BottomUpSubStructuredGraphSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::map< int, std::vector< PeriodicEntity > > PeriodicMap
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:172
std::shared_ptr< Geometry > GeometrySharedPtr
Definition: Geometry.h:53
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
std::size_t hash_range(Iter first, Iter last)
Definition: HashUtils.hpp:68
double NekDouble
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255