Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Assembly (e.g. local to global) base mapping routines
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
37 
38 using namespace std;
39 
40 namespace Nektar
41 {
42  namespace MultiRegions
43  {
44  /**
45  * @class AssemblyMap
46  * This class acts as a base class for constructing mappings between
47  * local, global and boundary degrees of freedom. It holds the storage
48  * for the maps and provides the accessors needed to retrieve them.
49  *
50  * There are two derived classes: AssemblyMapCG and
51  * AssemblyMapDG. These perform the actual construction of the
52  * maps within their specific contexts.
53  *
54  */
55 
56  /// Rounds a double precision number to an integer.
58  {
59  return int(x > 0.0 ? x + 0.5 : x - 0.5);
60  }
61 
62  /// Rounds an array of double precision numbers to integers.
64  {
65  int size = inarray.num_elements();
66  ASSERTL1(outarray.num_elements()>=size,"Array sizes not compatible");
67 
68  NekDouble x;
69  for(int i = 0; i < size; i++)
70  {
71  x = inarray[i];
72  outarray[i] = int(x > 0.0 ? x + 0.5 : x - 0.5);
73  }
74  }
75 
76  /**
77  * Initialises an empty mapping.
78  */
79  AssemblyMap::AssemblyMap():
80  m_session(),
81  m_comm(),
82  m_hash(0),
83  m_numLocalBndCoeffs(0),
84  m_numGlobalBndCoeffs(0),
85  m_numLocalDirBndCoeffs(0),
86  m_numGlobalDirBndCoeffs(0),
87  m_solnType(eNoSolnType),
88  m_bndSystemBandWidth(0),
89  m_successiveRHS(0),
90  m_gsh(0),
91  m_bndGsh(0)
92  {
93  }
94 
97  const std::string variable):
98  m_session(pSession),
99  m_comm(pSession->GetComm()),
100  m_hash(0),
101  m_numLocalBndCoeffs(0),
102  m_numGlobalBndCoeffs(0),
103  m_numLocalDirBndCoeffs(0),
104  m_numGlobalDirBndCoeffs(0),
105  m_bndSystemBandWidth(0),
106  m_successiveRHS(0),
107  m_gsh(0),
108  m_bndGsh(0)
109  {
110  // Default value from Solver Info
111  m_solnType = pSession->GetSolverInfoAsEnum<GlobalSysSolnType>(
112  "GlobalSysSoln");
113  m_preconType = pSession->GetSolverInfoAsEnum<PreconditionerType>(
114  "Preconditioner");
115 
116  // Override values with data from GlobalSysSolnInfo section
117  if(pSession->DefinesGlobalSysSolnInfo(variable, "GlobalSysSoln"))
118  {
119  std::string sysSoln = pSession->GetGlobalSysSolnInfo(variable,
120  "GlobalSysSoln");
121  m_solnType = pSession->GetValueAsEnum<GlobalSysSolnType>(
122  "GlobalSysSoln", sysSoln);
123  }
124 
125  if(pSession->DefinesGlobalSysSolnInfo(variable, "Preconditioner"))
126  {
127  std::string precon = pSession->GetGlobalSysSolnInfo(variable,
128  "Preconditioner");
129  m_preconType = pSession->GetValueAsEnum<PreconditionerType>(
130  "Preconditioner", precon);
131  }
132 
133  if(pSession->DefinesGlobalSysSolnInfo(variable,
134  "IterativeSolverTolerance"))
135  {
136  m_iterativeTolerance = boost::lexical_cast<NekDouble>(
137  pSession->GetGlobalSysSolnInfo(variable,
138  "IterativeSolverTolerance").c_str());
139  }
140  else
141  {
142  pSession->LoadParameter("IterativeSolverTolerance",
145  }
146 
147 
148  if(pSession->DefinesGlobalSysSolnInfo(variable,
149  "MaxIterations"))
150  {
151  m_maxIterations = boost::lexical_cast<int>(
152  pSession->GetGlobalSysSolnInfo(variable,
153  "MaxIterations").c_str());
154  }
155  else
156  {
157  pSession->LoadParameter("MaxIterations",
159  5000);
160  }
161 
162 
163  if(pSession->DefinesGlobalSysSolnInfo(variable,"SuccessiveRHS"))
164  {
165  m_successiveRHS = boost::lexical_cast<int>(
166  pSession->GetGlobalSysSolnInfo(variable,
167  "SuccessiveRHS").c_str());
168  }
169  else
170  {
171  pSession->LoadParameter("SuccessiveRHS",
172  m_successiveRHS,0);
173  }
174 
175  }
176 
177  /**
178  * Create a new level of mapping using the information in
179  * multiLevelGraph and performing the following steps:
180  */
182  AssemblyMap* oldLevelMap,
183  const BottomUpSubStructuredGraphSharedPtr& multiLevelGraph):
184  m_session(oldLevelMap->m_session),
185  m_comm(oldLevelMap->GetComm()),
186  m_hash(0),
187  m_solnType(oldLevelMap->m_solnType),
188  m_preconType(oldLevelMap->m_preconType),
189  m_maxIterations(oldLevelMap->m_maxIterations),
190  m_iterativeTolerance(oldLevelMap->m_iterativeTolerance),
191  m_successiveRHS(oldLevelMap->m_successiveRHS),
192  m_gsh(oldLevelMap->m_gsh),
193  m_bndGsh(oldLevelMap->m_bndGsh),
194  m_lowestStaticCondLevel(oldLevelMap->m_lowestStaticCondLevel)
195  {
196  int i;
197  int j;
198  int cnt;
199 
200  //--------------------------------------------------------------
201  // -- Extract information from the input argument
202  int numGlobalBndCoeffsOld = oldLevelMap->GetNumGlobalBndCoeffs();
203  int numGlobalDirBndCoeffsOld = oldLevelMap->GetNumGlobalDirBndCoeffs();
204  int numLocalBndCoeffsOld = oldLevelMap->GetNumLocalBndCoeffs();
205  int numLocalDirBndCoeffsOld = oldLevelMap->GetNumLocalDirBndCoeffs();
206  bool signChangeOld = oldLevelMap->GetSignChange();
207 
208  int staticCondLevelOld = oldLevelMap->GetStaticCondLevel();
209  int numPatchesOld = oldLevelMap->GetNumPatches();
210  GlobalSysSolnType solnTypeOld = oldLevelMap->GetGlobalSysSolnType();
211  const Array<OneD, const unsigned int>& numLocalBndCoeffsPerPatchOld = oldLevelMap->GetNumLocalBndCoeffsPerPatch();
212  //--------------------------------------------------------------
213 
214  //--------------------------------------------------------------
215  int newLevel = staticCondLevelOld+1;
216  /** - STEP 1: setup a mask array to determine to which patch
217  * of the new level every patch of the current
218  * level belongs. To do so we make four arrays,
219  * #gloPatchMask, #globHomPatchMask,
220  * #locPatchMask_NekDouble and #locPatchMask.
221  * These arrays are then used to check which local
222  * dofs of the old level belong to which patch of
223  * the new level
224  */
225  Array<OneD, NekDouble> globPatchMask (numGlobalBndCoeffsOld,-1.0);
226  Array<OneD, NekDouble> globHomPatchMask (globPatchMask+numGlobalDirBndCoeffsOld);
227  Array<OneD, NekDouble> locPatchMask_NekDouble(numLocalBndCoeffsOld,-3.0);
228  Array<OneD, int> locPatchMask (numLocalBndCoeffsOld);
229 
230  // Fill the array globPatchMask as follows:
231  // - The first part (i.e. the glob bnd dofs) is filled with the
232  // value -1
233  // - The second part (i.e. the glob interior dofs) is numbered
234  // according to the patch it belongs to (i.e. dofs in first block
235  // all are numbered 0, the second block numbered are 1, etc...)
236  multiLevelGraph->MaskPatches(newLevel,globHomPatchMask);
237 
238  // Map from Global Dofs to Local Dofs
239  // As a result, we know for each local dof whether
240  // it is mapped to the boundary of the next level, or to which
241  // patch. Based upon this, we can than later associate every patch
242  // of the current level with a patch in the next level.
243  oldLevelMap->GlobalToLocalBndWithoutSign(globPatchMask,locPatchMask_NekDouble);
244 
245  // Convert the result to an array of integers rather than doubles
246  RoundNekDoubleToInt(locPatchMask_NekDouble,locPatchMask);
247 
248  /** - STEP 2: We calculate how many local bnd dofs of the
249  * old level belong to the boundaries of each patch at
250  * the new level. We need this information to set up the
251  * mapping between different levels.
252  */
253 
254  // Retrieve the number of patches at the next level
255  int numPatchesWithIntNew = multiLevelGraph->GetNpatchesWithInterior(newLevel);
256  int numPatchesNew = numPatchesWithIntNew;
257 
258  // Allocate memory to store the number of local dofs associated to
259  // each of elemental boundaries of these patches
260  std::map<int, int> numLocalBndCoeffsPerPatchNew;
261  for(int i = 0; i < numPatchesNew; i++)
262  {
263  numLocalBndCoeffsPerPatchNew[i] = 0;
264  }
265 
266  int minval;
267  int maxval;
268  int curPatch;
269  for(i = cnt = 0; i < numPatchesOld; i++)
270  {
271  // For every patch at the current level, the mask array
272  // locPatchMask should be filled with
273  // - the same (positive) number for each entry (which will
274  // correspond to the patch at the next level it belongs to)
275  // - the same (positive) number for each entry, except some
276  // entries that are -1 (the enties correspond to -1, will be
277  // mapped to the local boundary of the next level patch given
278  // by the positive number)
279  // - -1 for all entries. In this case, we will make an
280  // additional patch only consisting of boundaries at the next
281  // level
282  minval = *min_element(&locPatchMask[cnt],
283  &locPatchMask[cnt]+numLocalBndCoeffsPerPatchOld[i]);
284  maxval = *max_element(&locPatchMask[cnt],
285  &locPatchMask[cnt]+numLocalBndCoeffsPerPatchOld[i]);
286  ASSERTL0((minval==maxval)||(minval==-1),"These values should never be the same");
287 
288  if(maxval == -1)
289  {
290  curPatch = numPatchesNew;
291  numLocalBndCoeffsPerPatchNew[curPatch] = 0;
292  numPatchesNew++;
293  }
294  else
295  {
296  curPatch = maxval;
297  }
298 
299  for(j = 0; j < numLocalBndCoeffsPerPatchOld[i]; j++ )
300  {
301  ASSERTL0((locPatchMask[cnt]==maxval)||(locPatchMask[cnt]==minval),
302  "These values should never be the same");
303  if(locPatchMask[cnt] == -1)
304  {
305  ++numLocalBndCoeffsPerPatchNew[curPatch];
306  }
307  cnt++;
308  }
309  }
310 
311  // Count how many local dofs of the old level are mapped
312  // to the local boundary dofs of the new level
314  m_numPatches = numLocalBndCoeffsPerPatchNew.size();
317  for(int i = 0; i < m_numPatches; i++)
318  {
319  m_numLocalBndCoeffsPerPatch[i] = (unsigned int) numLocalBndCoeffsPerPatchNew[i];
320  m_numLocalBndCoeffs += numLocalBndCoeffsPerPatchNew[i];
321  }
322  multiLevelGraph->GetNintDofsPerPatch(newLevel,m_numLocalIntCoeffsPerPatch);
323 
324  // Also initialise some more data members
325  m_solnType = solnTypeOld;
330  "This method should only be called for in "
331  "case of multi-level static condensation.");
332  m_staticCondLevel = newLevel;
333  m_signChange = signChangeOld;
334  m_numLocalDirBndCoeffs = numLocalDirBndCoeffsOld;
335  m_numGlobalDirBndCoeffs = numGlobalDirBndCoeffsOld;
336  m_numGlobalBndCoeffs = multiLevelGraph->GetInteriorOffset(newLevel) +
338  m_numGlobalCoeffs = multiLevelGraph->GetNumGlobalDofs(newLevel) +
341  if(m_signChange)
342  {
344  }
345 
347 
352 
353  // Set up an offset array that denotes the offset of the local
354  // boundary degrees of freedom of the next level
355  Array<OneD, int> numLocalBndCoeffsPerPatchOffset(m_numPatches+1,0);
356  for(int i = 1; i < m_numPatches+1; i++)
357  {
358  numLocalBndCoeffsPerPatchOffset[i] += numLocalBndCoeffsPerPatchOffset[i-1] + numLocalBndCoeffsPerPatchNew[i-1];
359  }
360 
361  int additionalPatchCnt = numPatchesWithIntNew;
362  int newid;
363  int blockid;
364  bool isBndDof;
365  NekDouble sign;
366  Array<OneD, int> bndDofPerPatchCnt(m_numPatches,0);
367  for(i = cnt = 0; i < numPatchesOld; i++)
368  {
369  minval = *min_element(&locPatchMask[cnt],&locPatchMask[cnt]+numLocalBndCoeffsPerPatchOld[i]);
370  maxval = *max_element(&locPatchMask[cnt],&locPatchMask[cnt]+numLocalBndCoeffsPerPatchOld[i]);
371  ASSERTL0((minval==maxval)||(minval==-1),"These values should never be the same");
372 
373  if(maxval == -1)
374  {
375  curPatch = additionalPatchCnt;
376  additionalPatchCnt++;
377  }
378  else
379  {
380  curPatch = maxval;
381  }
382 
383  for(j = 0; j < numLocalBndCoeffsPerPatchOld[i]; j++ )
384  {
385  ASSERTL0((locPatchMask[cnt]==maxval)||(locPatchMask[cnt]==minval),
386  "These values should never be the same");
387 
388  sign = oldLevelMap->GetLocalToGlobalBndSign(cnt);
389 
390  if(locPatchMask[cnt] == -1)
391  {
392  newid = numLocalBndCoeffsPerPatchOffset[curPatch];
393 
394  m_localToGlobalBndMap[newid] = oldLevelMap->GetLocalToGlobalBndMap(cnt);
395  if(m_signChange)
396  {
397  m_localToGlobalBndSign[ newid ] = sign;
398  }
399 
400  blockid = bndDofPerPatchCnt[curPatch];
401  isBndDof = true;
402 
403 
404  numLocalBndCoeffsPerPatchOffset[curPatch]++;
405  bndDofPerPatchCnt[curPatch]++;
406  }
407  else
408  {
409  newid = oldLevelMap->GetLocalToGlobalBndMap(cnt) -
410  m_numGlobalBndCoeffs+m_numLocalBndCoeffs;
411 
412  blockid = oldLevelMap->GetLocalToGlobalBndMap(cnt)-
413  m_numGlobalDirBndCoeffs - multiLevelGraph->GetInteriorOffset(newLevel,curPatch);
414  isBndDof = false;
415  }
416 
417  sign = isBndDof?1.0:sign;
418 
419  m_patchMapFromPrevLevel->SetPatchMap(cnt,curPatch,blockid,isBndDof,sign);
420 
421  cnt++;
422  }
423  }
425 
426  // Postprocess the computed information - Update the old
427  // level with the mapping to new evel
428  // oldLevelMap->SetLocalBndToNextLevelMap(oldLocalBndToNewLevelMap,oldLocalBndToNewLevelSign);
429  // - Construct the next level mapping object
430  if(m_staticCondLevel < (multiLevelGraph->GetNlevels()-1))
431  {
433  }
434  }
435 
436 
438  {
439  }
440 
441 
442  /**
443  * The bandwidth calculated corresponds to what is referred to as
444  * half-bandwidth. If the elements of the matrix are designated as
445  * a_ij, it corresponds to the maximum value of |i-j| for non-zero
446  * a_ij. As a result, the value also corresponds to the number of
447  * sub- or super-diagonals.
448  *
449  * The bandwith can be calculated elementally as it corresponds to the
450  * maximal elemental bandwith (i.e. the maximal difference in global
451  * DOF index for every element).
452  *
453  * We here calculate the bandwith of the global boundary system (as
454  * used for static condensation).
455  */
457  {
458  int i,j;
459  int cnt = 0;
460  int locSize;
461  int maxId;
462  int minId;
463  int bwidth = -1;
464  for(i = 0; i < m_numPatches; ++i)
465  {
466  locSize = m_numLocalBndCoeffsPerPatch[i];
467  maxId = -1;
468  minId = m_numLocalBndCoeffs+1;
469  for(j = 0; j < locSize; j++)
470  {
472  {
473  if(m_localToGlobalBndMap[cnt+j] > maxId)
474  {
475  maxId = m_localToGlobalBndMap[cnt+j];
476  }
477 
478  if(m_localToGlobalBndMap[cnt+j] < minId)
479  {
480  minId = m_localToGlobalBndMap[cnt+j];
481  }
482  }
483  }
484  bwidth = (bwidth>(maxId-minId))?bwidth:(maxId-minId);
485 
486  cnt+=locSize;
487  }
488 
489  m_bndSystemBandWidth = bwidth;
490  }
491 
492 
493  int AssemblyMap::v_GetLocalToGlobalMap(const int i) const
494  {
495  ASSERTL0(false, "Not defined for this type of mapping.");
496  return 0;
497  }
498 
500  {
501  ASSERTL0(false, "Not defined for this type of mapping.");
502  return 0;
503  }
504 
506  {
507  ASSERTL0(false, "Not defined for this type of mapping.");
508  return 0;
509  }
510 
512  {
513  ASSERTL0(false, "Not defined for this type of mapping.");
514  static Array<OneD,const int> result;
515  return result;
516  }
517 
519  {
520  ASSERTL0(false, "Not defined for this type of mapping.");
521  static Array<OneD, const int> result;
522  return result;
523  }
524 
526  {
527  ASSERTL0(false, "Not defined for this type of mapping.");
528  static Array<OneD, const int> result;
529  return result;
530  }
531 
533  {
534  ASSERTL0(false, "Not defined for this type of mapping.");
535  return 0.0;
536  }
537 
539  {
540  ASSERTL0(false, "Not defined for this type of mapping.");
541  static Array<OneD, NekDouble> result;
542  return result;
543  }
544 
546  const Array<OneD, const NekDouble>& loc,
547  Array<OneD, NekDouble>& global,
548  bool useComm) const
549  {
550  ASSERTL0(false, "Not defined for this type of mapping.");
551  }
552 
554  const NekVector<NekDouble>& loc,
555  NekVector< NekDouble>& global,
556  bool useComm) const
557  {
558  ASSERTL0(false, "Not defined for this type of mapping.");
559  }
560 
562  const Array<OneD, const NekDouble>& global,
563  Array<OneD, NekDouble>& loc) const
564  {
565  ASSERTL0(false, "Not defined for this type of mapping.");
566  }
567 
569  const NekVector<NekDouble>& global,
570  NekVector< NekDouble>& loc) const
571  {
572  ASSERTL0(false, "Not defined for this type of mapping.");
573  }
574 
576  const Array<OneD, const NekDouble> &loc,
577  Array<OneD, NekDouble> &global) const
578  {
579  ASSERTL0(false, "Not defined for this type of mapping.");
580  }
581 
583  const NekVector<NekDouble>& loc,
584  NekVector< NekDouble>& global) const
585  {
586  ASSERTL0(false, "Not defined for this type of mapping.");
587  }
588 
590  Array<OneD, NekDouble>& pGlobal) const
591  {
592  // Do nothing here since multi-level static condensation uses a
593  // AssemblyMap and thus will call this routine in serial.
594  }
595 
597  NekVector< NekDouble>& pGlobal) const
598  {
599  // Do nothing here since multi-level static condensation uses a
600  // AssemblyMap and thus will call this routine in serial.
601  }
602 
604  Array<OneD, NekDouble>& pGlobal,
605  int offset) const
606  {
607  // Do nothing here since multi-level static condensation uses a
608  // AssemblyMap and thus will call this routine in serial.
609  }
610 
612  {
613  ASSERTL0(false, "Not defined for this type of mapping.");
614  return 0;
615  }
616 
618  {
619  ASSERTL0(false, "Not defined for this type of mapping.");
620  return 0;
621  }
622 
624  {
625  ASSERTL0(false, "Not defined for this type of mapping.");
626  return 0;
627  }
628 
630  {
631  ASSERTL0(false, "Not defined for this type of mapping.");
632  return 0;
633  }
634 
636  {
637  ASSERTL0(false, "Not defined for this type of mapping.");
638  return 0;
639  }
640 
642  {
643  ASSERTL0(false, "Not defined for this type of mapping.");
644  return 0;
645  }
646 
648  {
649  ASSERTL0(false, "Not defined for this type of mapping.");
650  return 0;
651  }
652 
654  {
655  ASSERTL0(false, "Not defined for this type of mapping.");
656  return 0;
657  }
658 
660  {
661  ASSERTL0(false, "Not defined for this type of mapping.");
662  static Array<OneD, const int> result;
663  return result;
664  }
665 
666  boost::shared_ptr<AssemblyMap> AssemblyMap::v_LinearSpaceMap(
667  const ExpList &locexp, GlobalSysSolnType solnType)
668  {
669  ASSERTL0(false, "Not defined for this sub class");
670  static boost::shared_ptr<AssemblyMap> result;
671  return result;
672  }
673 
675  {
676  return m_comm;
677  }
678 
679  size_t AssemblyMap::GetHash() const
680  {
681  return m_hash;
682  }
683 
684  int AssemblyMap::GetLocalToGlobalMap(const int i) const
685  {
686  return v_GetLocalToGlobalMap(i);
687  }
688 
690  {
691  return v_GetGlobalToUniversalMap(i);
692  }
693 
695  {
697  }
698 
700  {
701  return v_GetLocalToGlobalMap();
702  }
703 
705  {
706  return v_GetGlobalToUniversalMap();
707  }
708 
710  {
712  }
713 
715  {
716  return v_GetLocalToGlobalSign(i);
717  }
718 
720  {
721  return v_GetLocalToGlobalSign();
722  }
723 
725  const Array<OneD, const NekDouble>& loc,
726  Array<OneD, NekDouble>& global,
727  bool useComm) const
728  {
729  v_LocalToGlobal(loc,global,useComm);
730  }
731 
733  const NekVector<NekDouble>& loc,
734  NekVector< NekDouble>& global,
735  bool useComm) const
736  {
737  v_LocalToGlobal(loc,global,useComm);
738  }
739 
741  const Array<OneD, const NekDouble>& global,
742  Array<OneD, NekDouble>& loc) const
743  {
744  v_GlobalToLocal(global,loc);
745  }
746 
748  const NekVector<NekDouble>& global,
749  NekVector< NekDouble>& loc) const
750  {
751  v_GlobalToLocal(global,loc);
752  }
753 
755  const Array<OneD, const NekDouble> &loc,
756  Array<OneD, NekDouble> &global) const
757  {
758  v_Assemble(loc,global);
759  }
760 
762  const NekVector<NekDouble>& loc,
763  NekVector< NekDouble>& global) const
764  {
765  v_Assemble(loc,global);
766  }
767 
769  Array<OneD, NekDouble>& pGlobal) const
770  {
771  v_UniversalAssemble(pGlobal);
772  }
773 
775  NekVector< NekDouble>& pGlobal) const
776  {
777  v_UniversalAssemble(pGlobal);
778  }
779 
781  Array<OneD, NekDouble>& pGlobal,
782  int offset) const
783  {
784  v_UniversalAssemble(pGlobal, offset);
785  }
786 
788  {
789  return v_GetFullSystemBandWidth();
790  }
791 
793  {
794  return v_GetNumNonDirVertexModes();
795  }
796 
798  {
799  return v_GetNumNonDirEdgeModes();
800  }
801 
803  {
804  return v_GetNumNonDirFaceModes();
805  }
806 
808  {
809  return v_GetNumDirEdges();
810  }
811 
813  {
814  return v_GetNumDirFaces();
815  }
816 
818  {
819  return v_GetNumNonDirEdges();
820  }
821 
823  {
824  return v_GetNumNonDirFaces();
825  }
826 
828  {
829  return v_GetExtraDirEdges();
830  }
831 
832  boost::shared_ptr<AssemblyMap> AssemblyMap::LinearSpaceMap(const ExpList &locexp, GlobalSysSolnType solnType)
833  {
834  return v_LinearSpaceMap(locexp, solnType);
835  }
836 
837  int AssemblyMap::GetLocalToGlobalBndMap(const int i) const
838  {
839  return m_localToGlobalBndMap[i];
840  }
841 
842  const Array<OneD,const int>&
844  {
845  return m_localToGlobalBndMap;
846  }
847 
849  {
850  return m_signChange;
851  }
852 
853 
856  {
857  return m_localToGlobalBndSign;
858  }
859 
861  {
863  }
864 
866  {
868  }
869 
871  {
872  if(m_signChange)
873  {
874  return m_localToGlobalBndSign[i];
875  }
876  else
877  {
878  return 1.0;
879  }
880  }
881 
882 
884  const int i)
885  {
887  }
888 
889 
891  const int i)
892  {
893  ASSERTL1(i < m_bndCondTraceToGlobalTraceMap.num_elements(),
894  "Index out of range.");
896  }
897 
900  {
901  return m_bndCondTraceToGlobalTraceMap;
902  }
903 
905  {
906  if(m_signChange)
907  {
909  }
910  else
911  {
912  return 1.0;
913  }
914  }
915 
916  const Array<OneD,const int>&
918  {
920  }
921 
922 
924  {
926  }
927 
928 
930  {
931  return m_numLocalDirBndCoeffs;
932  }
933 
935  {
936  return m_numLocalBndCoeffs;
937  }
938 
940  {
941  return m_numGlobalBndCoeffs;
942  }
943 
945  {
946  return m_numLocalCoeffs;
947  }
948 
950  {
951  return m_numGlobalCoeffs;
952  }
953 
955  {
956  return m_systemSingular;
957  }
958 
960  const NekVector<NekDouble>& global,
962  int offset) const
963  {
964  GlobalToLocalBnd(global.GetPtr(), loc.GetPtr(), offset);
965  }
966 
967 
969  const NekVector<NekDouble>& global,
970  NekVector<NekDouble>& loc) const
971  {
972  GlobalToLocalBnd(global.GetPtr(), loc.GetPtr());
973  }
974 
975 
977  const Array<OneD, const NekDouble>& global,
978  Array<OneD,NekDouble>& loc, int offset) const
979  {
980  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
981  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs-offset,"Global vector is not of correct dimension");
982 
983  // offset input data by length "offset" for Dirichlet boundary conditions.
985  Vmath::Vcopy(m_numGlobalBndCoeffs-offset, global.get(), 1, tmp.get() + offset, 1);
986 
987  if(m_signChange)
988  {
990  }
991  else
992  {
993  Vmath::Gathr(m_numLocalBndCoeffs, tmp.get(), m_localToGlobalBndMap.get(), loc.get());
994  }
995  }
996 
997 
999  const Array<OneD, const NekDouble>& global,
1000  Array<OneD,NekDouble>& loc) const
1001  {
1002  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1003  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs,"Global vector is not of correct dimension");
1004 
1005  if(m_signChange)
1006  {
1007  Vmath::Gathr(m_numLocalBndCoeffs, m_localToGlobalBndSign.get(), global.get(), m_localToGlobalBndMap.get(), loc.get());
1008  }
1009  else
1010  {
1011  Vmath::Gathr(m_numLocalBndCoeffs, global.get(), m_localToGlobalBndMap.get(), loc.get());
1012  }
1013  }
1014 
1016  const NekVector<NekDouble>& loc,
1017  NekVector<NekDouble>& global,
1018  int offset) const
1019  {
1020  LocalBndToGlobal(loc.GetPtr(), global.GetPtr(), offset);
1021  }
1022 
1023 
1025  const Array<OneD, const NekDouble>& loc,
1026  Array<OneD,NekDouble>& global,
1027  int offset) const
1028  {
1029  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1030  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs-offset,"Global vector is not of correct dimension");
1031 
1032  // offset input data by length "offset" for Dirichlet boundary conditions.
1034 
1035  if(m_signChange)
1036  {
1038  }
1039  else
1040  {
1041  Vmath::Scatr(m_numLocalBndCoeffs, loc.get(), m_localToGlobalBndMap.get(), tmp.get());
1042  }
1043 
1044  UniversalAssembleBnd(tmp);
1045  Vmath::Vcopy(m_numGlobalBndCoeffs-offset, tmp.get()+offset, 1, global.get(), 1);
1046  }
1047 
1049  const NekVector<NekDouble>& loc,
1050  NekVector<NekDouble>& global) const
1051  {
1052  LocalBndToGlobal(loc.GetPtr(), global.GetPtr());
1053  }
1054 
1055 
1057  const Array<OneD, const NekDouble>& loc,
1058  Array<OneD,NekDouble>& global) const
1059  {
1060  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1061  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs,"Global vector is not of correct dimension");
1062 
1063  if(m_signChange)
1064  {
1065  Vmath::Scatr(m_numLocalBndCoeffs, m_localToGlobalBndSign.get(), loc.get(), m_localToGlobalBndMap.get(), global.get());
1066  }
1067  else
1068  {
1069  Vmath::Scatr(m_numLocalBndCoeffs, loc.get(), m_localToGlobalBndMap.get(), global.get());
1070  }
1071  }
1072 
1073 
1075  const NekVector<NekDouble>& loc,
1076  NekVector<NekDouble>& global, int offset) const
1077  {
1078  AssembleBnd(loc.GetPtr(), global.GetPtr(), offset);
1079  }
1080 
1081 
1083  const NekVector<NekDouble>& loc,
1084  NekVector<NekDouble>& global) const
1085  {
1086  AssembleBnd(loc.GetPtr(), global.GetPtr());
1087  }
1088 
1089 
1091  const Array<OneD,const NekDouble>& loc,
1092  Array<OneD, NekDouble>& global, int offset) const
1093  {
1094  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local array is not of correct dimension");
1095  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs-offset,"Global array is not of correct dimension");
1097 
1098  if(m_signChange)
1099  {
1101  }
1102  else
1103  {
1104  Vmath::Assmb(m_numLocalBndCoeffs,loc.get(), m_localToGlobalBndMap.get(), tmp.get());
1105  }
1106  UniversalAssembleBnd(tmp);
1107  Vmath::Vcopy(m_numGlobalBndCoeffs-offset, tmp.get() + offset, 1, global.get(), 1);
1108  }
1109 
1110 
1112  const Array<OneD, const NekDouble>& loc,
1113  Array<OneD, NekDouble>& global) const
1114  {
1115  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1116  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs,"Global vector is not of correct dimension");
1117 
1118  Vmath::Zero(m_numGlobalBndCoeffs, global.get(), 1);
1119 
1120  if(m_signChange)
1121  {
1123  loc.get(), m_localToGlobalBndMap.get(), global.get());
1124  }
1125  else
1126  {
1127  Vmath::Assmb(m_numLocalBndCoeffs,loc.get(), m_localToGlobalBndMap.get(), global.get());
1128  }
1129  UniversalAssembleBnd(global);
1130  }
1131 
1133  Array<OneD, NekDouble>& pGlobal) const
1134  {
1135  ASSERTL1(pGlobal.num_elements() == m_numGlobalBndCoeffs,
1136  "Wrong size.");
1137  Gs::Gather(pGlobal, Gs::gs_add, m_bndGsh);
1138  }
1139 
1141  NekVector< NekDouble>& pGlobal) const
1142  {
1143  UniversalAssembleBnd(pGlobal.GetPtr());
1144  }
1145 
1147  Array<OneD, NekDouble>& pGlobal,
1148  int offset) const
1149  {
1150  Array<OneD, NekDouble> tmp(offset);
1151  if (offset > 0) Vmath::Vcopy(offset, pGlobal, 1, tmp, 1);
1152  UniversalAssembleBnd(pGlobal);
1153  if (offset > 0) Vmath::Vcopy(offset, tmp, 1, pGlobal, 1);
1154  }
1155 
1157  {
1158  return m_bndSystemBandWidth;
1159  }
1160 
1162  {
1163  return m_staticCondLevel;
1164  }
1165 
1167  {
1168  return m_numPatches;
1169  }
1170 
1173  {
1175  }
1176 
1177 
1180  {
1182  }
1183 
1184  const AssemblyMapSharedPtr
1186  {
1188  }
1189 
1190  const PatchMapSharedPtr&
1192  const
1193  {
1194  return m_patchMapFromPrevLevel;
1195  }
1196 
1198  {
1199  return !( m_nextLevelLocalToGlobalMap.get() );
1200  }
1201 
1202 
1204  {
1205  return m_solnType;
1206  }
1207 
1209  {
1210  return m_preconType;
1211  }
1212 
1214  {
1215  return m_iterativeTolerance;
1216  }
1217 
1219  {
1220  return m_maxIterations;
1221  }
1222 
1224  {
1225  return m_successiveRHS;
1226  }
1227 
1229  const Array<OneD, const NekDouble>& global,
1230  Array<OneD,NekDouble>& loc)
1231  {
1232  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1233  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs,"Global vector is not of correct dimension");
1234 
1235  Vmath::Gathr(m_numLocalBndCoeffs, global.get(), m_localToGlobalBndMap.get(), loc.get());
1236  }
1237 
1239  std::ostream &out, std::string variable, bool printHeader) const
1240  {
1242  = m_session->GetComm()->GetRowComm();
1243  bool isRoot = vRowComm->GetRank() == 0;
1244  int n = vRowComm->GetSize();
1245  int i;
1246 
1247  // Determine number of global degrees of freedom.
1248  int globBndCnt = 0, globDirCnt = 0;
1249 
1250  for (i = 0; i < m_numGlobalBndCoeffs; ++i)
1251  {
1253  {
1254  globBndCnt++;
1255 
1256  if (i < m_numGlobalDirBndCoeffs)
1257  {
1258  globDirCnt++;
1259  }
1260  }
1261  }
1262 
1263  int globCnt = m_numGlobalCoeffs - m_numGlobalBndCoeffs + globBndCnt;
1264 
1265  // Calculate maximum valency
1267  Array<OneD, NekDouble> tmpGlob(m_numGlobalBndCoeffs, 0.0);
1268 
1269  Vmath::Assmb(
1270  m_numLocalBndCoeffs, tmpLoc.get(),
1271  m_localToGlobalBndMap.get(), tmpGlob.get());
1272  UniversalAssembleBnd(tmpGlob);
1273 
1274  int totGlobDof = globCnt;
1275  int totGlobBndDof = globBndCnt;
1276  int totGlobDirDof = globDirCnt;
1277  int totLocalDof = m_numLocalCoeffs;
1278  int totLocalBndDof = m_numLocalBndCoeffs;
1279  int totLocalDirDof = m_numLocalDirBndCoeffs;
1280 
1281  int meanValence = 0;
1282  int maxValence = 0;
1283  int minValence = 10000000;
1284  for (int i = 0; i < m_numGlobalBndCoeffs; ++i)
1285  {
1287  {
1288  continue;
1289  }
1290 
1291  if (tmpGlob[i] > maxValence)
1292  {
1293  maxValence = tmpGlob[i];
1294  }
1295  if (tmpGlob[i] < minValence)
1296  {
1297  minValence = tmpGlob[i];
1298  }
1299  meanValence += tmpGlob[i];
1300  }
1301 
1302  vRowComm->AllReduce(maxValence, LibUtilities::ReduceMax);
1303  vRowComm->AllReduce(minValence, LibUtilities::ReduceMin);
1304  vRowComm->AllReduce(meanValence, LibUtilities::ReduceSum);
1305  vRowComm->AllReduce(totGlobDof, LibUtilities::ReduceSum);
1306  vRowComm->AllReduce(totGlobBndDof, LibUtilities::ReduceSum);
1307  vRowComm->AllReduce(totGlobDirDof, LibUtilities::ReduceSum);
1308  vRowComm->AllReduce(totLocalDof, LibUtilities::ReduceSum);
1309  vRowComm->AllReduce(totLocalBndDof, LibUtilities::ReduceSum);
1310  vRowComm->AllReduce(totLocalDirDof, LibUtilities::ReduceSum);
1311 
1312  meanValence /= totGlobBndDof;
1313 
1314  if (isRoot)
1315  {
1316  if (printHeader)
1317  {
1318  out << "Assembly map statistics for field " << variable
1319  << ":" << endl;
1320  }
1321 
1322  out << " - Number of local/global dof : "
1323  << totLocalDof << " " << totGlobDof << endl;
1324  out << " - Number of local/global boundary dof : "
1325  << totLocalBndDof << " " << totGlobBndDof << endl;
1326  out << " - Number of local/global Dirichlet dof : "
1327  << totLocalDirDof << " " << totGlobDirDof << endl;
1328  out << " - dof valency (min/max/mean) : "
1329  << minValence << " " << maxValence << " " << meanValence
1330  << endl;
1331 
1332  if (n > 1)
1333  {
1334  NekDouble mean = m_numLocalCoeffs, mean2 = mean * mean;
1335  NekDouble minval = mean, maxval = mean;
1336  Array<OneD, NekDouble> tmp(1);
1337 
1338  for (i = 1; i < n; ++i)
1339  {
1340  vRowComm->Recv(i, tmp);
1341  mean += tmp[0];
1342  mean2 += tmp[0]*tmp[0];
1343 
1344  if (tmp[0] > maxval)
1345  {
1346  maxval = tmp[0];
1347  }
1348  if (tmp[0] < minval)
1349  {
1350  minval = tmp[0];
1351  }
1352  }
1353 
1354  if (maxval > 0.1)
1355  {
1356  out << " - Local dof dist. (min/max/mean/dev) : "
1357  << minval << " " << maxval << " " << (mean / n)
1358  << " " << sqrt(mean2/n - mean*mean/n/n) << endl;
1359  }
1360 
1361  vRowComm->Block();
1362 
1363  mean = minval = maxval = m_numLocalBndCoeffs;
1364  mean2 = mean * mean;
1365 
1366  for (i = 1; i < n; ++i)
1367  {
1368  vRowComm->Recv(i, tmp);
1369  mean += tmp[0];
1370  mean2 += tmp[0]*tmp[0];
1371 
1372  if (tmp[0] > maxval)
1373  {
1374  maxval = tmp[0];
1375  }
1376  if (tmp[0] < minval)
1377  {
1378  minval = tmp[0];
1379  }
1380  }
1381 
1382  out << " - Local bnd dof dist. (min/max/mean/dev) : "
1383  << minval << " " << maxval << " " << (mean / n) << " "
1384  << sqrt(mean2/n - mean*mean/n/n) << endl;
1385  }
1386  }
1387  else
1388  {
1389  Array<OneD, NekDouble> tmp(1);
1390  tmp[0] = m_numLocalCoeffs;
1391  vRowComm->Send(0, tmp);
1392  vRowComm->Block();
1393  tmp[0] = m_numLocalBndCoeffs;
1394  vRowComm->Send(0, tmp);
1395  }
1396 
1397  // Either we have no more levels in the static condensation, or we
1398  // are not multi-level.
1400  {
1401  return;
1402  }
1403 
1404  int level = 2;
1406  while (tmp->m_nextLevelLocalToGlobalMap)
1407  {
1408  tmp = tmp->m_nextLevelLocalToGlobalMap;
1409  ++level;
1410  }
1411 
1412  // Print out multi-level static condensation information.
1413  if (n > 1)
1414  {
1415  if (isRoot)
1416  {
1417  NekDouble mean = level, mean2 = mean * mean;
1418  int minval = level, maxval = level;
1419 
1420  Array<OneD, NekDouble> tmpRecv(1);
1421  for (i = 1; i < n; ++i)
1422  {
1423  vRowComm->Recv(i, tmpRecv);
1424  mean += tmpRecv[0];
1425  mean2 += tmpRecv[0]*tmpRecv[0];
1426 
1427  if (tmpRecv[0] > maxval)
1428  {
1429  maxval = (int)(tmpRecv[0] + 0.5);
1430  }
1431  if (tmpRecv[0] < minval)
1432  {
1433  minval = (int)(tmpRecv[0] + 0.5);
1434  }
1435  }
1436 
1437  out << " - M-level sc. dist. (min/max/mean/dev) : "
1438  << minval << " " << maxval << " " << (mean / n) << " "
1439  << sqrt(mean2/n - mean*mean/n/n) << endl;
1440  }
1441  else
1442  {
1443  Array<OneD, NekDouble> tmpSend(1);
1444  tmpSend[0] = level;
1445  vRowComm->Send(0, tmpSend);
1446  }
1447  }
1448  else
1449  {
1450  out << " - Number of static cond. levels : "
1451  << level << endl;
1452  }
1453 
1454  if (isRoot)
1455  {
1456  out << "Stats at lowest static cond. level:" << endl;
1457  }
1458  tmp->PrintStats(out, variable, false);
1459  }
1460  } // namespace
1461 } // namespace
const Array< OneD, const int > & GetBndCondTraceToGlobalTraceMap()
void PrintStats(std::ostream &out, std::string variable, bool printHeader=true) const
PreconditionerType GetPreconType() const
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
bool m_systemSingular
Flag indicating if the system is singular or not.
Definition: AssemblyMap.h:322
const Array< OneD, const int > & GetGlobalToUniversalMap()
bool m_signChange
Flag indicating if modes require sign reversal.
Definition: AssemblyMap.h:347
int GetLocalToGlobalBndMap(const int i) const
Retrieve the global index of a given local boundary mode.
virtual int v_GetNumNonDirFaceModes() const
const Array< OneD, NekDouble > & GetLocalToGlobalSign() const
PreconditionerType m_preconType
Type type of preconditioner to use in iterative solver.
Definition: AssemblyMap.h:370
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:316
LibUtilities::CommSharedPtr GetComm()
Retrieves the communicator.
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: AssemblyMap.h:308
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:239
void Gathr(int n, const T *x, const int *y, T *z)
Gather vector z[i] = x[y[i]].
Definition: Vmath.cpp:644
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:27
virtual void v_LocalToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, bool useComm) const
boost::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:53
int GetNumLocalCoeffs() const
Returns the total number of local coefficients.
virtual int v_GetNumNonDirVertexModes() const
static const NekDouble kNekIterativeTol
Array< OneD, int > m_bndCondTraceToGlobalTraceMap
Integer map of bnd cond trace number to global trace number.
Definition: AssemblyMap.h:358
int GetBndSystemBandWidth() const
Returns the bandwidth of the boundary system.
virtual int v_GetNumNonDirFaces() const
size_t GetHash() const
Retrieves the hash of this map.
boost::shared_ptr< BottomUpSubStructuredGraph > BottomUpSubStructuredGraphSharedPtr
STL namespace.
bool AtLastLevel() const
Returns true if this is the last level in the multi-level static condensation.
int GetNumPatches() const
Returns the number of patches in this static condensation level.
Array< OneD, const NekDouble > GetLocalToGlobalBndSign() const
Retrieve the sign change for all local boundary modes.
NekDouble GetIterativeTolerance() const
int m_numLocalCoeffs
Total number of local coefficients.
Definition: AssemblyMap.h:333
const Array< OneD, const int > & GetExtraDirEdges()
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
PatchMapSharedPtr m_patchMapFromPrevLevel
Mapping information for previous level in MultiLevel Solver.
Definition: AssemblyMap.h:410
virtual void v_GlobalToLocal(const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const
boost::shared_ptr< PatchMap > PatchMapSharedPtr
const Array< OneD, const int > & GetGlobalToUniversalBndMapUnique()
virtual void v_UniversalAssemble(Array< OneD, NekDouble > &pGlobal) const
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:397
Base class for constructing local to global mapping of degrees of freedom.
Definition: AssemblyMap.h:59
int m_successiveRHS
sucessive RHS for iterative solver
Definition: AssemblyMap.h:379
size_t m_hash
Hash for map.
Definition: AssemblyMap.h:311
void UniversalAssemble(Array< OneD, NekDouble > &pGlobal) const
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:101
boost::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
void LocalToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, bool useComm=true) const
int GetNumGlobalDirBndCoeffs() const
Returns the number of global Dirichlet boundary coefficients.
virtual boost::shared_ptr< AssemblyMap > v_LinearSpaceMap(const ExpList &locexp, GlobalSysSolnType solnType)
Generate a linear space mapping from existing mapping.
GlobalSysSolnType GetGlobalSysSolnType() const
Returns the method of solving global systems.
void RoundNekDoubleToInt(const Array< OneD, const NekDouble > inarray, Array< OneD, int > outarray)
Rounds an array of double precision numbers to integers.
Definition: AssemblyMap.cpp:63
int m_bndSystemBandWidth
The bandwith of the global bnd system.
Definition: AssemblyMap.h:367
NekDouble m_iterativeTolerance
Tolerance for iterative solver.
Definition: AssemblyMap.h:376
Array< OneD, unsigned int > m_numLocalBndCoeffsPerPatch
The number of bnd dofs per patch.
Definition: AssemblyMap.h:390
boost::shared_ptr< AssemblyMap > LinearSpaceMap(const ExpList &locexp, GlobalSysSolnType solnType)
GlobalSysSolnType m_solnType
The solution type of the global system.
Definition: AssemblyMap.h:365
Array< OneD, NekDouble > m_bndCondCoeffsToGlobalCoeffsSign
Integer map of bnd cond coeffs to global coefficients.
Definition: AssemblyMap.h:356
int m_numGlobalDirBndCoeffs
Number of Global Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:320
NekDouble GetBndCondCoeffsToGlobalCoeffsSign(const int i)
Returns the modal sign associated with a given boundary expansion mode.
int GetStaticCondLevel() const
Returns the level of static condensation for this map.
void GlobalToLocal(const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const
virtual int v_GetNumDirFaces() const
virtual int v_GetFullSystemBandWidth() const
void Scatr(int n, const T *x, const int *y, T *z)
Scatter vector z[y[i]] = x[i].
Definition: Vmath.cpp:673
const PatchMapSharedPtr & GetPatchMapFromPrevLevel(void) const
Returns the patch map from the previous level of the multi-level static condensation.
double NekDouble
void Assmb(int n, const T *x, const int *y, T *z)
Assemble z[y[i]] += x[i]; z should be zero'd first.
Definition: Vmath.cpp:709
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
The number of int dofs per patch.
Definition: AssemblyMap.h:392
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:350
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:318
bool GetSignChange()
Returns true if using a modal expansion requiring a change of sign of some modes. ...
const Array< OneD, const int > & GetLocalToGlobalMap()
virtual void v_Assemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
const Array< OneD, const int > & GetGlobalToUniversalBndMap()
virtual int v_GetNumNonDirEdgeModes() const
Array< OneD, int > m_bndCondCoeffsToGlobalCoeffsMap
Integer map of bnd cond coeffs to global coefficients.
Definition: AssemblyMap.h:354
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:314
int m_staticCondLevel
The level of recursion in the case of multi-level static condensation.
Definition: AssemblyMap.h:386
Array< OneD, NekDouble > m_localToGlobalBndSign
Integer sign of local boundary coeffs to global space.
Definition: AssemblyMap.h:352
int GetNumGlobalBndCoeffs() const
Returns the total number of global boundary coefficients.
AssemblyMap()
Default constructor.
Definition: AssemblyMap.cpp:79
int m_maxIterations
Maximum iterations for iterative solver.
Definition: AssemblyMap.h:373
const Array< OneD, const unsigned int > & GetNumLocalBndCoeffsPerPatch()
Returns the number of local boundary coefficients in each patch.
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:57
const Array< OneD, const int > & GetGlobalToUniversalMapUnique()
NekDouble GetLocalToGlobalBndSign(const int i) const
Retrieve the sign change of a given local boundary mode.
int GetNumLocalDirBndCoeffs() const
Returns the number of local Dirichlet boundary coefficients.
LibUtilities::SessionReaderSharedPtr m_session
Session object.
Definition: AssemblyMap.h:305
void LocalBndToGlobal(const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const
void AssembleBnd(const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const
virtual const Array< OneD, NekDouble > & v_GetLocalToGlobalSign() const
virtual const Array< OneD, const int > & v_GetExtraDirEdges()
int GetNumGlobalCoeffs() const
Returns the total number of global coefficients.
int GetNumLocalBndCoeffs() const
Returns the total number of local boundary coefficients.
Array< OneD, int > m_globalToUniversalBndMap
Integer map of process coeffs to universal space.
Definition: AssemblyMap.h:360
Array< OneD, int > m_globalToUniversalBndMapUnique
Integer map of unique process coeffs to universal space (signed)
Definition: AssemblyMap.h:362
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMapUnique()
No Solution type specified.
void UniversalAssembleBnd(Array< OneD, NekDouble > &pGlobal) const
int m_numGlobalCoeffs
Total number of global coefficients.
Definition: AssemblyMap.h:344
void GlobalToLocalBnd(const NekVector< NekDouble > &global, NekVector< NekDouble > &loc, int offset) const
const Array< OneD, const unsigned int > & GetNumLocalIntCoeffsPerPatch()
Returns the number of local interior coefficients in each patch.
virtual int v_GetNumNonDirEdges() const
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:373
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:228
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
const AssemblyMapSharedPtr GetNextLevelLocalToGlobalMap() const
Returns the local to global mapping for the next level in the multi-level static condensation.
void Assemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
Array< OneD, DataType > & GetPtr()
Definition: NekVector.cpp:230
virtual const Array< OneD, const int > & v_GetLocalToGlobalMap()
virtual int v_GetNumDirEdges() const
int m_numPatches
The number of patches (~elements) in the current level.
Definition: AssemblyMap.h:388
bool GetSingularSystem() const
Retrieves if the system is singular (true) or not (false)