Nektar++
AssemblyMap.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File AssemblyMap.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: Assembly (e.g. local to global) base mapping routines
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <boost/core/ignore_unused.hpp>
36 
38 
39 using namespace std;
40 
41 namespace Nektar
42 {
43  namespace MultiRegions
44  {
45  /**
46  * @class AssemblyMap
47  * This class acts as a base class for constructing mappings between
48  * local, global and boundary degrees of freedom. It holds the storage
49  * for the maps and provides the accessors needed to retrieve them.
50  *
51  * There are two derived classes: AssemblyMapCG and
52  * AssemblyMapDG. These perform the actual construction of the
53  * maps within their specific contexts.
54  *
55  */
56 
57  /// Rounds a double precision number to an integer.
59  {
60  return int(x > 0.0 ? x + 0.5 : x - 0.5);
61  }
62 
63  /// Rounds an array of double precision numbers to integers.
65  {
66  int size = inarray.num_elements();
67  ASSERTL1(outarray.num_elements()>=size,"Array sizes not compatible");
68 
69  NekDouble x;
70  for(int i = 0; i < size; i++)
71  {
72  x = inarray[i];
73  outarray[i] = int(x > 0.0 ? x + 0.5 : x - 0.5);
74  }
75  }
76 
77  /**
78  * Initialises an empty mapping.
79  */
80  AssemblyMap::AssemblyMap():
81  m_session(),
82  m_comm(),
83  m_hash(0),
84  m_numLocalBndCoeffs(0),
85  m_numGlobalBndCoeffs(0),
86  m_numLocalDirBndCoeffs(0),
87  m_numGlobalDirBndCoeffs(0),
88  m_solnType(eNoSolnType),
89  m_bndSystemBandWidth(0),
90  m_successiveRHS(0),
91  m_gsh(0),
92  m_bndGsh(0)
93  {
94  }
95 
98  const std::string variable):
99  m_session(pSession),
100  m_comm(pSession->GetComm()),
101  m_hash(0),
107  m_successiveRHS(0),
108  m_gsh(0),
109  m_bndGsh(0)
110  {
111  // Default value from Solver Info
112  m_solnType = pSession->GetSolverInfoAsEnum<GlobalSysSolnType>(
113  "GlobalSysSoln");
114  m_preconType = pSession->GetSolverInfoAsEnum<PreconditionerType>(
115  "Preconditioner");
116 
117  // Override values with data from GlobalSysSolnInfo section
118  if(pSession->DefinesGlobalSysSolnInfo(variable, "GlobalSysSoln"))
119  {
120  std::string sysSoln = pSession->GetGlobalSysSolnInfo(variable,
121  "GlobalSysSoln");
122  m_solnType = pSession->GetValueAsEnum<GlobalSysSolnType>(
123  "GlobalSysSoln", sysSoln);
124  }
125 
126  if(pSession->DefinesGlobalSysSolnInfo(variable, "Preconditioner"))
127  {
128  std::string precon = pSession->GetGlobalSysSolnInfo(variable,
129  "Preconditioner");
130  m_preconType = pSession->GetValueAsEnum<PreconditionerType>(
131  "Preconditioner", precon);
132  }
133 
134  if(pSession->DefinesGlobalSysSolnInfo(variable,
135  "IterativeSolverTolerance"))
136  {
137  m_iterativeTolerance = boost::lexical_cast<NekDouble>(
138  pSession->GetGlobalSysSolnInfo(variable,
139  "IterativeSolverTolerance").c_str());
140  }
141  else
142  {
143  pSession->LoadParameter("IterativeSolverTolerance",
146  }
147 
148 
149  if(pSession->DefinesGlobalSysSolnInfo(variable,
150  "MaxIterations"))
151  {
152  m_maxIterations = boost::lexical_cast<int>(
153  pSession->GetGlobalSysSolnInfo(variable,
154  "MaxIterations").c_str());
155  }
156  else
157  {
158  pSession->LoadParameter("MaxIterations",
160  5000);
161  }
162 
163 
164  if(pSession->DefinesGlobalSysSolnInfo(variable,"SuccessiveRHS"))
165  {
166  m_successiveRHS = boost::lexical_cast<int>(
167  pSession->GetGlobalSysSolnInfo(variable,
168  "SuccessiveRHS").c_str());
169  }
170  else
171  {
172  pSession->LoadParameter("SuccessiveRHS",
173  m_successiveRHS,0);
174  }
175 
176  }
177 
178  /**
179  * Create a new level of mapping using the information in
180  * multiLevelGraph and performing the following steps:
181  */
183  AssemblyMap* oldLevelMap,
184  const BottomUpSubStructuredGraphSharedPtr& multiLevelGraph):
185  m_session(oldLevelMap->m_session),
186  m_comm(oldLevelMap->GetComm()),
187  m_hash(0),
188  m_solnType(oldLevelMap->m_solnType),
189  m_preconType(oldLevelMap->m_preconType),
190  m_maxIterations(oldLevelMap->m_maxIterations),
192  m_successiveRHS(oldLevelMap->m_successiveRHS),
193  m_gsh(oldLevelMap->m_gsh),
194  m_bndGsh(oldLevelMap->m_bndGsh),
196  {
197  int i;
198  int j;
199  int cnt;
200 
201  //--------------------------------------------------------------
202  // -- Extract information from the input argument
203  int numGlobalBndCoeffsOld = oldLevelMap->GetNumGlobalBndCoeffs();
204  int numGlobalDirBndCoeffsOld = oldLevelMap->GetNumGlobalDirBndCoeffs();
205  int numLocalBndCoeffsOld = oldLevelMap->GetNumLocalBndCoeffs();
206  int numLocalDirBndCoeffsOld = oldLevelMap->GetNumLocalDirBndCoeffs();
207  bool signChangeOld = oldLevelMap->GetSignChange();
208 
209  int staticCondLevelOld = oldLevelMap->GetStaticCondLevel();
210  int numPatchesOld = oldLevelMap->GetNumPatches();
211  GlobalSysSolnType solnTypeOld = oldLevelMap->GetGlobalSysSolnType();
212  const Array<OneD, const unsigned int>& numLocalBndCoeffsPerPatchOld = oldLevelMap->GetNumLocalBndCoeffsPerPatch();
213  //--------------------------------------------------------------
214 
215  //--------------------------------------------------------------
216  int newLevel = staticCondLevelOld+1;
217  /** - STEP 1: setup a mask array to determine to which patch
218  * of the new level every patch of the current
219  * level belongs. To do so we make four arrays,
220  * #gloPatchMask, #globHomPatchMask,
221  * #locPatchMask_NekDouble and #locPatchMask.
222  * These arrays are then used to check which local
223  * dofs of the old level belong to which patch of
224  * the new level
225  */
226  Array<OneD, NekDouble> globPatchMask (numGlobalBndCoeffsOld,-1.0);
227  Array<OneD, NekDouble> globHomPatchMask (globPatchMask+numGlobalDirBndCoeffsOld);
228  Array<OneD, NekDouble> locPatchMask_NekDouble(numLocalBndCoeffsOld,-3.0);
229  Array<OneD, int> locPatchMask (numLocalBndCoeffsOld);
230 
231  // Fill the array globPatchMask as follows:
232  // - The first part (i.e. the glob bnd dofs) is filled with the
233  // value -1
234  // - The second part (i.e. the glob interior dofs) is numbered
235  // according to the patch it belongs to (i.e. dofs in first block
236  // all are numbered 0, the second block numbered are 1, etc...)
237  multiLevelGraph->MaskPatches(newLevel,globHomPatchMask);
238 
239  // Map from Global Dofs to Local Dofs
240  // As a result, we know for each local dof whether
241  // it is mapped to the boundary of the next level, or to which
242  // patch. Based upon this, we can than later associate every patch
243  // of the current level with a patch in the next level.
244  oldLevelMap->GlobalToLocalBndWithoutSign(globPatchMask,locPatchMask_NekDouble);
245 
246  // Convert the result to an array of integers rather than doubles
247  RoundNekDoubleToInt(locPatchMask_NekDouble,locPatchMask);
248 
249  /** - STEP 2: We calculate how many local bnd dofs of the
250  * old level belong to the boundaries of each patch at
251  * the new level. We need this information to set up the
252  * mapping between different levels.
253  */
254 
255  // Retrieve the number of patches at the next level
256  int numPatchesWithIntNew = multiLevelGraph->GetNpatchesWithInterior(newLevel);
257  int numPatchesNew = numPatchesWithIntNew;
258 
259  // Allocate memory to store the number of local dofs associated to
260  // each of elemental boundaries of these patches
261  std::map<int, int> numLocalBndCoeffsPerPatchNew;
262  for(int i = 0; i < numPatchesNew; i++)
263  {
264  numLocalBndCoeffsPerPatchNew[i] = 0;
265  }
266 
267  int minval;
268  int maxval;
269  int curPatch;
270  for(i = cnt = 0; i < numPatchesOld; i++)
271  {
272  // For every patch at the current level, the mask array
273  // locPatchMask should be filled with
274  // - the same (positive) number for each entry (which will
275  // correspond to the patch at the next level it belongs to)
276  // - the same (positive) number for each entry, except some
277  // entries that are -1 (the enties correspond to -1, will be
278  // mapped to the local boundary of the next level patch given
279  // by the positive number)
280  // - -1 for all entries. In this case, we will make an
281  // additional patch only consisting of boundaries at the next
282  // level
283  minval = *min_element(&locPatchMask[cnt],
284  &locPatchMask[cnt]+numLocalBndCoeffsPerPatchOld[i]);
285  maxval = *max_element(&locPatchMask[cnt],
286  &locPatchMask[cnt]+numLocalBndCoeffsPerPatchOld[i]);
287  ASSERTL0((minval==maxval)||(minval==-1),"These values should never be the same");
288 
289  if(maxval == -1)
290  {
291  curPatch = numPatchesNew;
292  numLocalBndCoeffsPerPatchNew[curPatch] = 0;
293  numPatchesNew++;
294  }
295  else
296  {
297  curPatch = maxval;
298  }
299 
300  for(j = 0; j < numLocalBndCoeffsPerPatchOld[i]; j++ )
301  {
302  ASSERTL0((locPatchMask[cnt]==maxval)||(locPatchMask[cnt]==minval),
303  "These values should never be the same");
304  if(locPatchMask[cnt] == -1)
305  {
306  ++numLocalBndCoeffsPerPatchNew[curPatch];
307  }
308  cnt++;
309  }
310  }
311 
312  // Count how many local dofs of the old level are mapped
313  // to the local boundary dofs of the new level
315  m_numPatches = numLocalBndCoeffsPerPatchNew.size();
318  for(int i = 0; i < m_numPatches; i++)
319  {
320  m_numLocalBndCoeffsPerPatch[i] = (unsigned int) numLocalBndCoeffsPerPatchNew[i];
321  m_numLocalBndCoeffs += numLocalBndCoeffsPerPatchNew[i];
322  }
323  multiLevelGraph->GetNintDofsPerPatch(newLevel,m_numLocalIntCoeffsPerPatch);
324 
325  // Also initialise some more data members
326  m_solnType = solnTypeOld;
331  "This method should only be called for in "
332  "case of multi-level static condensation.");
333  m_staticCondLevel = newLevel;
334  m_signChange = signChangeOld;
335  m_numLocalDirBndCoeffs = numLocalDirBndCoeffsOld;
336  m_numGlobalDirBndCoeffs = numGlobalDirBndCoeffsOld;
337  m_numGlobalBndCoeffs = multiLevelGraph->GetInteriorOffset(newLevel) +
339  m_numGlobalCoeffs = multiLevelGraph->GetNumGlobalDofs(newLevel) +
342  if(m_signChange)
343  {
345  }
346 
348 
353 
354  // Set up an offset array that denotes the offset of the local
355  // boundary degrees of freedom of the next level
356  Array<OneD, int> numLocalBndCoeffsPerPatchOffset(m_numPatches+1,0);
357  for(int i = 1; i < m_numPatches+1; i++)
358  {
359  numLocalBndCoeffsPerPatchOffset[i] += numLocalBndCoeffsPerPatchOffset[i-1] + numLocalBndCoeffsPerPatchNew[i-1];
360  }
361 
362  int additionalPatchCnt = numPatchesWithIntNew;
363  int newid;
364  int blockid;
365  bool isBndDof;
366  NekDouble sign;
367  Array<OneD, int> bndDofPerPatchCnt(m_numPatches,0);
368  for(i = cnt = 0; i < numPatchesOld; i++)
369  {
370  minval = *min_element(&locPatchMask[cnt],&locPatchMask[cnt]+numLocalBndCoeffsPerPatchOld[i]);
371  maxval = *max_element(&locPatchMask[cnt],&locPatchMask[cnt]+numLocalBndCoeffsPerPatchOld[i]);
372  ASSERTL0((minval==maxval)||(minval==-1),"These values should never be the same");
373 
374  if(maxval == -1)
375  {
376  curPatch = additionalPatchCnt;
377  additionalPatchCnt++;
378  }
379  else
380  {
381  curPatch = maxval;
382  }
383 
384  for(j = 0; j < numLocalBndCoeffsPerPatchOld[i]; j++ )
385  {
386  ASSERTL0((locPatchMask[cnt]==maxval)||(locPatchMask[cnt]==minval),
387  "These values should never be the same");
388 
389  sign = oldLevelMap->GetLocalToGlobalBndSign(cnt);
390 
391  if(locPatchMask[cnt] == -1)
392  {
393  newid = numLocalBndCoeffsPerPatchOffset[curPatch];
394 
395  m_localToGlobalBndMap[newid] = oldLevelMap->GetLocalToGlobalBndMap(cnt);
396  if(m_signChange)
397  {
398  m_localToGlobalBndSign[ newid ] = sign;
399  }
400 
401  blockid = bndDofPerPatchCnt[curPatch];
402  isBndDof = true;
403 
404 
405  numLocalBndCoeffsPerPatchOffset[curPatch]++;
406  bndDofPerPatchCnt[curPatch]++;
407  }
408  else
409  {
410  newid = oldLevelMap->GetLocalToGlobalBndMap(cnt) -
411  m_numGlobalBndCoeffs+m_numLocalBndCoeffs;
412 
413  blockid = oldLevelMap->GetLocalToGlobalBndMap(cnt)-
414  m_numGlobalDirBndCoeffs - multiLevelGraph->GetInteriorOffset(newLevel,curPatch);
415  isBndDof = false;
416  }
417 
418  sign = isBndDof?1.0:sign;
419 
420  m_patchMapFromPrevLevel->SetPatchMap(cnt,curPatch,blockid,isBndDof,sign);
421 
422  cnt++;
423  }
424  }
426 
427  // Postprocess the computed information - Update the old
428  // level with the mapping to new evel
429  // oldLevelMap->SetLocalBndToNextLevelMap(oldLocalBndToNewLevelMap,oldLocalBndToNewLevelSign);
430  // - Construct the next level mapping object
431  if(m_staticCondLevel < (multiLevelGraph->GetNlevels()-1))
432  {
434  }
435  }
436 
437 
439  {
440  }
441 
442 
443  /**
444  * The bandwidth calculated corresponds to what is referred to as
445  * half-bandwidth. If the elements of the matrix are designated as
446  * a_ij, it corresponds to the maximum value of |i-j| for non-zero
447  * a_ij. As a result, the value also corresponds to the number of
448  * sub- or super-diagonals.
449  *
450  * The bandwith can be calculated elementally as it corresponds to the
451  * maximal elemental bandwith (i.e. the maximal difference in global
452  * DOF index for every element).
453  *
454  * We here calculate the bandwith of the global boundary system (as
455  * used for static condensation).
456  */
458  {
459  int i,j;
460  int cnt = 0;
461  int locSize;
462  int maxId;
463  int minId;
464  int bwidth = -1;
465  for(i = 0; i < m_numPatches; ++i)
466  {
467  locSize = m_numLocalBndCoeffsPerPatch[i];
468  maxId = -1;
469  minId = m_numLocalBndCoeffs+1;
470  for(j = 0; j < locSize; j++)
471  {
473  {
474  if(m_localToGlobalBndMap[cnt+j] > maxId)
475  {
476  maxId = m_localToGlobalBndMap[cnt+j];
477  }
478 
479  if(m_localToGlobalBndMap[cnt+j] < minId)
480  {
481  minId = m_localToGlobalBndMap[cnt+j];
482  }
483  }
484  }
485  bwidth = (bwidth>(maxId-minId))?bwidth:(maxId-minId);
486 
487  cnt+=locSize;
488  }
489 
490  m_bndSystemBandWidth = bwidth;
491  }
492 
493 
494  int AssemblyMap::v_GetLocalToGlobalMap(const int i) const
495  {
496  boost::ignore_unused(i);
498  "Not defined for this type of mapping.");
499  return 0;
500  }
501 
503  {
504  boost::ignore_unused(i);
506  "Not defined for this type of mapping.");
507  return 0;
508  }
509 
511  {
512  boost::ignore_unused(i);
514  "Not defined for this type of mapping.");
515  return 0;
516  }
517 
519  {
521  "Not defined for this type of mapping.");
522  static Array<OneD,const int> result;
523  return result;
524  }
525 
527  {
529  "Not defined for this type of mapping.");
530  static Array<OneD, const int> result;
531  return result;
532  }
533 
535  {
537  "Not defined for this type of mapping.");
538  static Array<OneD, const int> result;
539  return result;
540  }
541 
543  {
544  boost::ignore_unused(i);
546  "Not defined for this type of mapping.");
547  return 0.0;
548  }
549 
551  {
553  "Not defined for this type of mapping.");
554  static Array<OneD, NekDouble> result;
555  return result;
556  }
557 
560  Array<OneD, NekDouble>& global,
561  bool useComm) const
562  {
563  boost::ignore_unused(loc, global, useComm);
565  "Not defined for this type of mapping.");
566  }
567 
569  const NekVector<NekDouble>& loc,
570  NekVector< NekDouble>& global,
571  bool useComm) const
572  {
573  boost::ignore_unused(loc, global, useComm);
575  "Not defined for this type of mapping.");
576  }
577 
579  const Array<OneD, const NekDouble>& global,
581  {
582  boost::ignore_unused(loc, global);
584  "Not defined for this type of mapping.");
585  }
586 
588  const NekVector<NekDouble>& global,
589  NekVector< NekDouble>& loc) const
590  {
591  boost::ignore_unused(loc, global);
593  "Not defined for this type of mapping.");
594  }
595 
598  Array<OneD, NekDouble> &global) const
599  {
600  boost::ignore_unused(loc, global);
602  "Not defined for this type of mapping.");
603  }
604 
606  const NekVector<NekDouble>& loc,
607  NekVector< NekDouble>& global) const
608  {
609  boost::ignore_unused(loc, global);
611  "Not defined for this type of mapping.");
612  }
613 
615  Array<OneD, NekDouble>& pGlobal) const
616  {
617  boost::ignore_unused(pGlobal);
618  // Do nothing here since multi-level static condensation uses a
619  // AssemblyMap and thus will call this routine in serial.
620  }
621 
623  NekVector< NekDouble>& pGlobal) const
624  {
625  boost::ignore_unused(pGlobal);
626  // Do nothing here since multi-level static condensation uses a
627  // AssemblyMap and thus will call this routine in serial.
628  }
629 
631  Array<OneD, NekDouble>& pGlobal,
632  int offset) const
633  {
634  boost::ignore_unused(pGlobal, offset);
635  // Do nothing here since multi-level static condensation uses a
636  // AssemblyMap and thus will call this routine in serial.
637  }
638 
640  {
642  "Not defined for this type of mapping.");
643  return 0;
644  }
645 
647  {
649  "Not defined for this type of mapping.");
650  return 0;
651  }
652 
654  {
656  "Not defined for this type of mapping.");
657  return 0;
658  }
659 
661  {
663  "Not defined for this type of mapping.");
664  return 0;
665  }
666 
668  {
670  "Not defined for this type of mapping.");
671  return 0;
672  }
673 
675  {
677  "Not defined for this type of mapping.");
678  return 0;
679  }
680 
682  {
684  "Not defined for this type of mapping.");
685  return 0;
686  }
687 
689  {
691  "Not defined for this type of mapping.");
692  return 0;
693  }
694 
696  {
698  "Not defined for this type of mapping.");
699  static Array<OneD, const int> result;
700  return result;
701  }
702 
703  std::shared_ptr<AssemblyMap> AssemblyMap::v_LinearSpaceMap(
704  const ExpList &locexp, GlobalSysSolnType solnType)
705  {
706  boost::ignore_unused(locexp, solnType);
708  "Not defined for this type of mapping.");
709  static std::shared_ptr<AssemblyMap> result;
710  return result;
711  }
712 
714  {
715  return m_comm;
716  }
717 
718  size_t AssemblyMap::GetHash() const
719  {
720  return m_hash;
721  }
722 
723  int AssemblyMap::GetLocalToGlobalMap(const int i) const
724  {
725  return v_GetLocalToGlobalMap(i);
726  }
727 
729  {
730  return v_GetGlobalToUniversalMap(i);
731  }
732 
734  {
736  }
737 
739  {
740  return v_GetLocalToGlobalMap();
741  }
742 
744  {
745  return v_GetGlobalToUniversalMap();
746  }
747 
749  {
751  }
752 
754  {
755  return v_GetLocalToGlobalSign(i);
756  }
757 
759  {
760  return v_GetLocalToGlobalSign();
761  }
762 
765  Array<OneD, NekDouble>& global,
766  bool useComm) const
767  {
768  v_LocalToGlobal(loc,global,useComm);
769  }
770 
772  const NekVector<NekDouble>& loc,
773  NekVector< NekDouble>& global,
774  bool useComm) const
775  {
776  v_LocalToGlobal(loc,global,useComm);
777  }
778 
780  const Array<OneD, const NekDouble>& global,
782  {
783  v_GlobalToLocal(global,loc);
784  }
785 
787  const NekVector<NekDouble>& global,
788  NekVector< NekDouble>& loc) const
789  {
790  v_GlobalToLocal(global,loc);
791  }
792 
795  Array<OneD, NekDouble> &global) const
796  {
797  v_Assemble(loc,global);
798  }
799 
801  const NekVector<NekDouble>& loc,
802  NekVector< NekDouble>& global) const
803  {
804  v_Assemble(loc,global);
805  }
806 
808  Array<OneD, NekDouble>& pGlobal) const
809  {
810  v_UniversalAssemble(pGlobal);
811  }
812 
814  NekVector< NekDouble>& pGlobal) const
815  {
816  v_UniversalAssemble(pGlobal);
817  }
818 
820  Array<OneD, NekDouble>& pGlobal,
821  int offset) const
822  {
823  v_UniversalAssemble(pGlobal, offset);
824  }
825 
827  {
828  return v_GetFullSystemBandWidth();
829  }
830 
832  {
833  return v_GetNumNonDirVertexModes();
834  }
835 
837  {
838  return v_GetNumNonDirEdgeModes();
839  }
840 
842  {
843  return v_GetNumNonDirFaceModes();
844  }
845 
847  {
848  return v_GetNumDirEdges();
849  }
850 
852  {
853  return v_GetNumDirFaces();
854  }
855 
857  {
858  return v_GetNumNonDirEdges();
859  }
860 
862  {
863  return v_GetNumNonDirFaces();
864  }
865 
867  {
868  return v_GetExtraDirEdges();
869  }
870 
871  std::shared_ptr<AssemblyMap> AssemblyMap::LinearSpaceMap(const ExpList &locexp, GlobalSysSolnType solnType)
872  {
873  return v_LinearSpaceMap(locexp, solnType);
874  }
875 
876  int AssemblyMap::GetLocalToGlobalBndMap(const int i) const
877  {
878  return m_localToGlobalBndMap[i];
879  }
880 
881  const Array<OneD,const int>&
883  {
884  return m_localToGlobalBndMap;
885  }
886 
888  {
889  return m_signChange;
890  }
891 
892 
895  {
896  return m_localToGlobalBndSign;
897  }
898 
900  {
902  }
903 
905  {
907  }
908 
910  {
911  if(m_signChange)
912  {
913  return m_localToGlobalBndSign[i];
914  }
915  else
916  {
917  return 1.0;
918  }
919  }
920 
921 
923  const int i)
924  {
926  }
927 
928 
930  const int i)
931  {
932  ASSERTL1(i < m_bndCondTraceToGlobalTraceMap.num_elements(),
933  "Index out of range.");
935  }
936 
939  {
941  }
942 
944  {
945  if(m_signChange)
946  {
948  }
949  else
950  {
951  return 1.0;
952  }
953  }
954 
955  const Array<OneD,const int>&
957  {
959  }
960 
961 
963  {
965  }
966 
967 
969  {
970  return m_numLocalDirBndCoeffs;
971  }
972 
974  {
975  return m_numLocalBndCoeffs;
976  }
977 
979  {
980  return m_numGlobalBndCoeffs;
981  }
982 
984  {
985  return m_numLocalCoeffs;
986  }
987 
989  {
990  return m_numGlobalCoeffs;
991  }
992 
994  {
995  return m_systemSingular;
996  }
997 
999  const NekVector<NekDouble>& global,
1001  int offset) const
1002  {
1003  GlobalToLocalBnd(global.GetPtr(), loc.GetPtr(), offset);
1004  }
1005 
1006 
1008  const NekVector<NekDouble>& global,
1009  NekVector<NekDouble>& loc) const
1010  {
1011  GlobalToLocalBnd(global.GetPtr(), loc.GetPtr());
1012  }
1013 
1014 
1016  const Array<OneD, const NekDouble>& global,
1017  Array<OneD,NekDouble>& loc, int offset) const
1018  {
1019  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1020  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs-offset,"Global vector is not of correct dimension");
1021 
1022  // offset input data by length "offset" for Dirichlet boundary conditions.
1024  Vmath::Vcopy(m_numGlobalBndCoeffs-offset, global.get(), 1, tmp.get() + offset, 1);
1025 
1026  if(m_signChange)
1027  {
1029  }
1030  else
1031  {
1032  Vmath::Gathr(m_numLocalBndCoeffs, tmp.get(), m_localToGlobalBndMap.get(), loc.get());
1033  }
1034  }
1035 
1036 
1038  const Array<OneD, const NekDouble>& global,
1039  Array<OneD,NekDouble>& loc) const
1040  {
1041  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1042  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs,"Global vector is not of correct dimension");
1043 
1044  if(m_signChange)
1045  {
1046  Vmath::Gathr(m_numLocalBndCoeffs, m_localToGlobalBndSign.get(), global.get(), m_localToGlobalBndMap.get(), loc.get());
1047  }
1048  else
1049  {
1050  Vmath::Gathr(m_numLocalBndCoeffs, global.get(), m_localToGlobalBndMap.get(), loc.get());
1051  }
1052  }
1053 
1055  const NekVector<NekDouble>& loc,
1056  NekVector<NekDouble>& global,
1057  int offset) const
1058  {
1059  LocalBndToGlobal(loc.GetPtr(), global.GetPtr(), offset);
1060  }
1061 
1062 
1065  Array<OneD,NekDouble>& global,
1066  int offset) const
1067  {
1068  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1069  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs-offset,"Global vector is not of correct dimension");
1070 
1071  // offset input data by length "offset" for Dirichlet boundary conditions.
1073 
1074  if(m_signChange)
1075  {
1077  }
1078  else
1079  {
1080  Vmath::Scatr(m_numLocalBndCoeffs, loc.get(), m_localToGlobalBndMap.get(), tmp.get());
1081  }
1082 
1083  UniversalAssembleBnd(tmp);
1084  Vmath::Vcopy(m_numGlobalBndCoeffs-offset, tmp.get()+offset, 1, global.get(), 1);
1085  }
1086 
1088  const NekVector<NekDouble>& loc,
1089  NekVector<NekDouble>& global) const
1090  {
1091  LocalBndToGlobal(loc.GetPtr(), global.GetPtr());
1092  }
1093 
1094 
1097  Array<OneD,NekDouble>& global) const
1098  {
1099  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1100  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs,"Global vector is not of correct dimension");
1101 
1102  if(m_signChange)
1103  {
1104  Vmath::Scatr(m_numLocalBndCoeffs, m_localToGlobalBndSign.get(), loc.get(), m_localToGlobalBndMap.get(), global.get());
1105  }
1106  else
1107  {
1108  Vmath::Scatr(m_numLocalBndCoeffs, loc.get(), m_localToGlobalBndMap.get(), global.get());
1109  }
1110  }
1111 
1112 
1114  const NekVector<NekDouble>& loc,
1115  NekVector<NekDouble>& global, int offset) const
1116  {
1117  AssembleBnd(loc.GetPtr(), global.GetPtr(), offset);
1118  }
1119 
1120 
1122  const NekVector<NekDouble>& loc,
1123  NekVector<NekDouble>& global) const
1124  {
1125  AssembleBnd(loc.GetPtr(), global.GetPtr());
1126  }
1127 
1128 
1131  Array<OneD, NekDouble>& global, int offset) const
1132  {
1133  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local array is not of correct dimension");
1134  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs-offset,"Global array is not of correct dimension");
1136 
1137  if(m_signChange)
1138  {
1140  }
1141  else
1142  {
1143  Vmath::Assmb(m_numLocalBndCoeffs,loc.get(), m_localToGlobalBndMap.get(), tmp.get());
1144  }
1145  UniversalAssembleBnd(tmp);
1146  Vmath::Vcopy(m_numGlobalBndCoeffs-offset, tmp.get() + offset, 1, global.get(), 1);
1147  }
1148 
1149 
1152  Array<OneD, NekDouble>& global) const
1153  {
1154  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1155  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs,"Global vector is not of correct dimension");
1156 
1157  Vmath::Zero(m_numGlobalBndCoeffs, global.get(), 1);
1158 
1159  if(m_signChange)
1160  {
1162  loc.get(), m_localToGlobalBndMap.get(), global.get());
1163  }
1164  else
1165  {
1166  Vmath::Assmb(m_numLocalBndCoeffs,loc.get(), m_localToGlobalBndMap.get(), global.get());
1167  }
1168  UniversalAssembleBnd(global);
1169  }
1170 
1172  Array<OneD, NekDouble>& pGlobal) const
1173  {
1174  ASSERTL1(pGlobal.num_elements() == m_numGlobalBndCoeffs,
1175  "Wrong size.");
1176  Gs::Gather(pGlobal, Gs::gs_add, m_bndGsh);
1177  }
1178 
1180  NekVector< NekDouble>& pGlobal) const
1181  {
1182  UniversalAssembleBnd(pGlobal.GetPtr());
1183  }
1184 
1186  Array<OneD, NekDouble>& pGlobal,
1187  int offset) const
1188  {
1189  Array<OneD, NekDouble> tmp(offset);
1190  if (offset > 0) Vmath::Vcopy(offset, pGlobal, 1, tmp, 1);
1191  UniversalAssembleBnd(pGlobal);
1192  if (offset > 0) Vmath::Vcopy(offset, tmp, 1, pGlobal, 1);
1193  }
1194 
1196  {
1197  return m_bndSystemBandWidth;
1198  }
1199 
1201  {
1202  return m_staticCondLevel;
1203  }
1204 
1206  {
1207  return m_numPatches;
1208  }
1209 
1212  {
1214  }
1215 
1216 
1219  {
1221  }
1222 
1223  const AssemblyMapSharedPtr
1225  {
1227  }
1228 
1229  const PatchMapSharedPtr&
1231  const
1232  {
1233  return m_patchMapFromPrevLevel;
1234  }
1235 
1237  {
1238  return !( m_nextLevelLocalToGlobalMap.get() );
1239  }
1240 
1241 
1243  {
1244  return m_solnType;
1245  }
1246 
1248  {
1249  return m_preconType;
1250  }
1251 
1253  {
1254  return m_iterativeTolerance;
1255  }
1256 
1258  {
1259  return m_maxIterations;
1260  }
1261 
1263  {
1264  return m_successiveRHS;
1265  }
1266 
1268  const Array<OneD, const NekDouble>& global,
1270  {
1271  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1272  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs,"Global vector is not of correct dimension");
1273 
1274  Vmath::Gathr(m_numLocalBndCoeffs, global.get(), m_localToGlobalBndMap.get(), loc.get());
1275  }
1276 
1278  std::ostream &out, std::string variable, bool printHeader) const
1279  {
1281  = m_session->GetComm()->GetRowComm();
1282  bool isRoot = vRowComm->GetRank() == 0;
1283  int n = vRowComm->GetSize();
1284  int i;
1285 
1286  // Determine number of global degrees of freedom.
1287  int globBndCnt = 0, globDirCnt = 0;
1288 
1289  for (i = 0; i < m_numGlobalBndCoeffs; ++i)
1290  {
1292  {
1293  globBndCnt++;
1294 
1295  if (i < m_numGlobalDirBndCoeffs)
1296  {
1297  globDirCnt++;
1298  }
1299  }
1300  }
1301 
1302  int globCnt = m_numGlobalCoeffs - m_numGlobalBndCoeffs + globBndCnt;
1303 
1304  // Calculate maximum valency
1306  Array<OneD, NekDouble> tmpGlob(m_numGlobalBndCoeffs, 0.0);
1307 
1308  Vmath::Assmb(
1309  m_numLocalBndCoeffs, tmpLoc.get(),
1310  m_localToGlobalBndMap.get(), tmpGlob.get());
1311  UniversalAssembleBnd(tmpGlob);
1312 
1313  int totGlobDof = globCnt;
1314  int totGlobBndDof = globBndCnt;
1315  int totGlobDirDof = globDirCnt;
1316  int totLocalDof = m_numLocalCoeffs;
1317  int totLocalBndDof = m_numLocalBndCoeffs;
1318  int totLocalDirDof = m_numLocalDirBndCoeffs;
1319 
1320  int meanValence = 0;
1321  int maxValence = 0;
1322  int minValence = 10000000;
1323  for (int i = 0; i < m_numGlobalBndCoeffs; ++i)
1324  {
1326  {
1327  continue;
1328  }
1329 
1330  if (tmpGlob[i] > maxValence)
1331  {
1332  maxValence = tmpGlob[i];
1333  }
1334  if (tmpGlob[i] < minValence)
1335  {
1336  minValence = tmpGlob[i];
1337  }
1338  meanValence += tmpGlob[i];
1339  }
1340 
1341  vRowComm->AllReduce(maxValence, LibUtilities::ReduceMax);
1342  vRowComm->AllReduce(minValence, LibUtilities::ReduceMin);
1343  vRowComm->AllReduce(meanValence, LibUtilities::ReduceSum);
1344  vRowComm->AllReduce(totGlobDof, LibUtilities::ReduceSum);
1345  vRowComm->AllReduce(totGlobBndDof, LibUtilities::ReduceSum);
1346  vRowComm->AllReduce(totGlobDirDof, LibUtilities::ReduceSum);
1347  vRowComm->AllReduce(totLocalDof, LibUtilities::ReduceSum);
1348  vRowComm->AllReduce(totLocalBndDof, LibUtilities::ReduceSum);
1349  vRowComm->AllReduce(totLocalDirDof, LibUtilities::ReduceSum);
1350 
1351  meanValence /= totGlobBndDof;
1352 
1353  if (isRoot)
1354  {
1355  if (printHeader)
1356  {
1357  out << "Assembly map statistics for field " << variable
1358  << ":" << endl;
1359  }
1360 
1361  out << " - Number of local/global dof : "
1362  << totLocalDof << " " << totGlobDof << endl;
1363  out << " - Number of local/global boundary dof : "
1364  << totLocalBndDof << " " << totGlobBndDof << endl;
1365  out << " - Number of local/global Dirichlet dof : "
1366  << totLocalDirDof << " " << totGlobDirDof << endl;
1367  out << " - dof valency (min/max/mean) : "
1368  << minValence << " " << maxValence << " " << meanValence
1369  << endl;
1370 
1371  if (n > 1)
1372  {
1373  NekDouble mean = m_numLocalCoeffs, mean2 = mean * mean;
1374  NekDouble minval = mean, maxval = mean;
1375  Array<OneD, NekDouble> tmp(1);
1376 
1377  for (i = 1; i < n; ++i)
1378  {
1379  vRowComm->Recv(i, tmp);
1380  mean += tmp[0];
1381  mean2 += tmp[0]*tmp[0];
1382 
1383  if (tmp[0] > maxval)
1384  {
1385  maxval = tmp[0];
1386  }
1387  if (tmp[0] < minval)
1388  {
1389  minval = tmp[0];
1390  }
1391  }
1392 
1393  if (maxval > 0.1)
1394  {
1395  out << " - Local dof dist. (min/max/mean/dev) : "
1396  << minval << " " << maxval << " " << (mean / n)
1397  << " " << sqrt(mean2/n - mean*mean/n/n) << endl;
1398  }
1399 
1400  vRowComm->Block();
1401 
1402  mean = minval = maxval = m_numLocalBndCoeffs;
1403  mean2 = mean * mean;
1404 
1405  for (i = 1; i < n; ++i)
1406  {
1407  vRowComm->Recv(i, tmp);
1408  mean += tmp[0];
1409  mean2 += tmp[0]*tmp[0];
1410 
1411  if (tmp[0] > maxval)
1412  {
1413  maxval = tmp[0];
1414  }
1415  if (tmp[0] < minval)
1416  {
1417  minval = tmp[0];
1418  }
1419  }
1420 
1421  out << " - Local bnd dof dist. (min/max/mean/dev) : "
1422  << minval << " " << maxval << " " << (mean / n) << " "
1423  << sqrt(mean2/n - mean*mean/n/n) << endl;
1424  }
1425  }
1426  else
1427  {
1428  Array<OneD, NekDouble> tmp(1);
1429  tmp[0] = m_numLocalCoeffs;
1430  vRowComm->Send(0, tmp);
1431  vRowComm->Block();
1432  tmp[0] = m_numLocalBndCoeffs;
1433  vRowComm->Send(0, tmp);
1434  }
1435 
1436  // Either we have no more levels in the static condensation, or we
1437  // are not multi-level.
1439  {
1440  return;
1441  }
1442 
1443  int level = 2;
1445  while (tmp->m_nextLevelLocalToGlobalMap)
1446  {
1447  tmp = tmp->m_nextLevelLocalToGlobalMap;
1448  ++level;
1449  }
1450 
1451  // Print out multi-level static condensation information.
1452  if (n > 1)
1453  {
1454  if (isRoot)
1455  {
1456  NekDouble mean = level, mean2 = mean * mean;
1457  int minval = level, maxval = level;
1458 
1459  Array<OneD, NekDouble> tmpRecv(1);
1460  for (i = 1; i < n; ++i)
1461  {
1462  vRowComm->Recv(i, tmpRecv);
1463  mean += tmpRecv[0];
1464  mean2 += tmpRecv[0]*tmpRecv[0];
1465 
1466  if (tmpRecv[0] > maxval)
1467  {
1468  maxval = (int)(tmpRecv[0] + 0.5);
1469  }
1470  if (tmpRecv[0] < minval)
1471  {
1472  minval = (int)(tmpRecv[0] + 0.5);
1473  }
1474  }
1475 
1476  out << " - M-level sc. dist. (min/max/mean/dev) : "
1477  << minval << " " << maxval << " " << (mean / n) << " "
1478  << sqrt(mean2/n - mean*mean/n/n) << endl;
1479  }
1480  else
1481  {
1482  Array<OneD, NekDouble> tmpSend(1);
1483  tmpSend[0] = level;
1484  vRowComm->Send(0, tmpSend);
1485  }
1486  }
1487  else
1488  {
1489  out << " - Number of static cond. levels : "
1490  << level << endl;
1491  }
1492 
1493  if (isRoot)
1494  {
1495  out << "Stats at lowest static cond. level:" << endl;
1496  }
1497  tmp->PrintStats(out, variable, false);
1498  }
1499  } // namespace
1500 } // namespace
const Array< OneD, const int > & GetBndCondTraceToGlobalTraceMap()
void GlobalToLocalBnd(const NekVector< NekDouble > &global, NekVector< NekDouble > &loc, int offset) const
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
int GetNumLocalDirBndCoeffs() const
Returns the number of local Dirichlet boundary coefficients.
virtual int v_GetFullSystemBandWidth() const
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:209
bool m_systemSingular
Flag indicating if the system is singular or not.
Definition: AssemblyMap.h:321
const Array< OneD, const int > & GetGlobalToUniversalMap()
bool m_signChange
Flag indicating if modes require sign reversal.
Definition: AssemblyMap.h:346
PreconditionerType GetPreconType() const
PreconditionerType m_preconType
Type type of preconditioner to use in iterative solver.
Definition: AssemblyMap.h:369
NekDouble GetLocalToGlobalBndSign(const int i) const
Retrieve the sign change of a given local boundary mode.
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:315
int GetNumGlobalCoeffs() const
Returns the total number of global coefficients.
LibUtilities::CommSharedPtr GetComm()
Retrieves the communicator.
bool GetSingularSystem() const
Retrieves if the system is singular (true) or not (false)
int GetBndSystemBandWidth() const
Returns the bandwidth of the boundary system.
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: AssemblyMap.h:307
virtual ~AssemblyMap()
Destructor.
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
void Gathr(int n, const T *x, const int *y, T *z)
Gather vector z[i] = x[y[i]].
Definition: Vmath.cpp:647
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:16
static const NekDouble kNekIterativeTol
Array< OneD, int > m_bndCondTraceToGlobalTraceMap
Integer map of bnd cond trace number to global trace number.
Definition: AssemblyMap.h:357
virtual int v_GetNumNonDirFaceModes() const
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:53
virtual int v_GetNumNonDirEdgeModes() const
virtual int v_GetNumNonDirEdges() const
virtual void v_LocalToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, bool useComm) const
STL namespace.
int m_numLocalCoeffs
Total number of local coefficients.
Definition: AssemblyMap.h:332
const Array< OneD, const int > & GetExtraDirEdges()
void UniversalAssemble(Array< OneD, NekDouble > &pGlobal) const
PatchMapSharedPtr m_patchMapFromPrevLevel
Mapping information for previous level in MultiLevel Solver.
Definition: AssemblyMap.h:409
const Array< OneD, const int > & GetGlobalToUniversalBndMapUnique()
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:52
const Array< OneD, const int > & GetLocalToGlobalBndMap()
Retrieve the global indices of the local boundary modes.
AssemblyMapSharedPtr m_nextLevelLocalToGlobalMap
Map from the patches of the previous level to the patches of the current level.
Definition: AssemblyMap.h:396
virtual std::shared_ptr< AssemblyMap > v_LinearSpaceMap(const ExpList &locexp, GlobalSysSolnType solnType)
Generate a linear space mapping from existing mapping.
Base class for constructing local to global mapping of degrees of freedom.
Definition: AssemblyMap.h:58
virtual int v_GetNumDirEdges() const
int m_successiveRHS
sucessive RHS for iterative solver
Definition: AssemblyMap.h:378
const PatchMapSharedPtr & GetPatchMapFromPrevLevel(void) const
Returns the patch map from the previous level of the multi-level static condensation.
virtual void v_UniversalAssemble(Array< OneD, NekDouble > &pGlobal) const
int GetNumLocalBndCoeffs() const
Returns the total number of local boundary coefficients.
bool AtLastLevel() const
Returns true if this is the last level in the multi-level static condensation.
size_t m_hash
Hash for map.
Definition: AssemblyMap.h:310
int GetLocalToGlobalBndMap(const int i) const
Retrieve the global index of a given local boundary mode.
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:103
std::shared_ptr< PatchMap > PatchMapSharedPtr
int GetNumGlobalBndCoeffs() const
Returns the total number of global boundary coefficients.
void RoundNekDoubleToInt(const Array< OneD, const NekDouble > inarray, Array< OneD, int > outarray)
Rounds an array of double precision numbers to integers.
Definition: AssemblyMap.cpp:64
virtual void v_Assemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
int m_bndSystemBandWidth
The bandwith of the global bnd system.
Definition: AssemblyMap.h:366
virtual void v_GlobalToLocal(const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const
NekDouble m_iterativeTolerance
Tolerance for iterative solver.
Definition: AssemblyMap.h:375
virtual int v_GetNumNonDirFaces() const
NekDouble GetIterativeTolerance() 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
Array< OneD, NekDouble > m_bndCondCoeffsToGlobalCoeffsSign
Integer map of bnd cond coeffs to global coefficients.
Definition: AssemblyMap.h:355
int m_numGlobalDirBndCoeffs
Number of Global Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:319
NekDouble GetBndCondCoeffsToGlobalCoeffsSign(const int i)
Returns the modal sign associated with a given boundary expansion mode.
void LocalBndToGlobal(const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void Scatr(int n, const T *x, const int *y, T *z)
Scatter vector z[y[i]] = x[i].
Definition: Vmath.cpp:676
Array< OneD, const NekDouble > GetLocalToGlobalBndSign() const
Retrieve the sign change for all local boundary modes.
void UniversalAssembleBnd(Array< OneD, NekDouble > &pGlobal) const
size_t GetHash() const
Retrieves the hash of this map.
virtual int v_GetNumDirFaces() const
double NekDouble
void Assmb(int n, const T *x, const int *y, T *z)
Assemble z[y[i]] += x[i]; z should be zero&#39;d first.
Definition: Vmath.cpp:712
GlobalSysSolnType GetGlobalSysSolnType() const
Returns the method of solving global systems.
virtual int v_GetNumNonDirVertexModes() const
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
The number of int dofs per patch.
Definition: AssemblyMap.h:391
const AssemblyMapSharedPtr GetNextLevelLocalToGlobalMap() const
Returns the local to global mapping for the next level in the multi-level static condensation.
int GetNumGlobalDirBndCoeffs() const
Returns the number of global Dirichlet boundary coefficients.
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
const Array< OneD, const int > & GetBndCondCoeffsToGlobalCoeffsMap()
Retrieves the global indices corresponding to the boundary expansion modes.
int m_numLocalDirBndCoeffs
Number of Local Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:317
bool GetSignChange()
Returns true if using a modal expansion requiring a change of sign of some modes. ...
const Array< OneD, const int > & GetLocalToGlobalMap()
const Array< OneD, const int > & GetGlobalToUniversalBndMap()
Array< OneD, int > m_bndCondCoeffsToGlobalCoeffsMap
Integer map of bnd cond coeffs to global coefficients.
Definition: AssemblyMap.h:353
const Array< OneD, NekDouble > & GetLocalToGlobalSign() const
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:313
void LocalToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, bool useComm=true) const
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
void Assemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
AssemblyMap()
Default constructor.
Definition: AssemblyMap.cpp:80
int m_maxIterations
Maximum iterations for iterative solver.
Definition: AssemblyMap.h:372
int GetStaticCondLevel() const
Returns the level of static condensation for this map.
virtual const Array< OneD, NekDouble > & v_GetLocalToGlobalSign() const
const Array< OneD, const unsigned int > & GetNumLocalBndCoeffsPerPatch()
Returns the number of local boundary coefficients in each patch.
void AssembleBnd(const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const
void GlobalToLocalBndWithoutSign(const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc)
int RoundNekDoubleToInt(NekDouble x)
Rounds a double precision number to an integer.
Definition: AssemblyMap.cpp:58
const Array< OneD, const int > & GetGlobalToUniversalMapUnique()
int GetNumPatches() const
Returns the number of patches in this static condensation level.
void PrintStats(std::ostream &out, std::string variable, bool printHeader=true) const
LibUtilities::SessionReaderSharedPtr m_session
Session object.
Definition: AssemblyMap.h:304
int GetNumLocalCoeffs() const
Returns the total number of local coefficients.
virtual const Array< OneD, const int > & v_GetExtraDirEdges()
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
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMapUnique()
No Solution type specified.
int m_numGlobalCoeffs
Total number of global coefficients.
Definition: AssemblyMap.h:343
const Array< OneD, const unsigned int > & GetNumLocalIntCoeffsPerPatch()
Returns the number of local interior coefficients in each patch.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:376
std::shared_ptr< AssemblyMap > LinearSpaceMap(const ExpList &locexp, GlobalSysSolnType solnType)
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMap()
#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
Array< OneD, DataType > & GetPtr()
Definition: NekVector.cpp:227
virtual const Array< OneD, const int > & v_GetLocalToGlobalMap()
int m_numPatches
The number of patches (~elements) in the current level.
Definition: AssemblyMap.h:387
void GlobalToLocal(const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const