Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator 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) const
548  {
549  ASSERTL0(false, "Not defined for this type of mapping.");
550  }
551 
553  const NekVector<NekDouble>& loc,
554  NekVector< NekDouble>& global) const
555  {
556  ASSERTL0(false, "Not defined for this type of mapping.");
557  }
558 
560  const Array<OneD, const NekDouble>& global,
561  Array<OneD, NekDouble>& loc) const
562  {
563  ASSERTL0(false, "Not defined for this type of mapping.");
564  }
565 
567  const NekVector<NekDouble>& global,
568  NekVector< NekDouble>& loc) const
569  {
570  ASSERTL0(false, "Not defined for this type of mapping.");
571  }
572 
574  const Array<OneD, const NekDouble> &loc,
575  Array<OneD, NekDouble> &global) const
576  {
577  ASSERTL0(false, "Not defined for this type of mapping.");
578  }
579 
581  const NekVector<NekDouble>& loc,
582  NekVector< NekDouble>& global) const
583  {
584  ASSERTL0(false, "Not defined for this type of mapping.");
585  }
586 
588  Array<OneD, NekDouble>& pGlobal) const
589  {
590  // Do nothing here since multi-level static condensation uses a
591  // AssemblyMap and thus will call this routine in serial.
592  }
593 
595  NekVector< NekDouble>& pGlobal) const
596  {
597  // Do nothing here since multi-level static condensation uses a
598  // AssemblyMap and thus will call this routine in serial.
599  }
600 
602  Array<OneD, NekDouble>& pGlobal,
603  int offset) const
604  {
605  // Do nothing here since multi-level static condensation uses a
606  // AssemblyMap and thus will call this routine in serial.
607  }
608 
610  {
611  ASSERTL0(false, "Not defined for this type of mapping.");
612  return 0;
613  }
614 
616  {
617  ASSERTL0(false, "Not defined for this type of mapping.");
618  return 0;
619  }
620 
622  {
623  ASSERTL0(false, "Not defined for this type of mapping.");
624  return 0;
625  }
626 
628  {
629  ASSERTL0(false, "Not defined for this type of mapping.");
630  return 0;
631  }
632 
634  {
635  ASSERTL0(false, "Not defined for this type of mapping.");
636  return 0;
637  }
638 
640  {
641  ASSERTL0(false, "Not defined for this type of mapping.");
642  return 0;
643  }
644 
646  {
647  ASSERTL0(false, "Not defined for this type of mapping.");
648  return 0;
649  }
650 
652  {
653  ASSERTL0(false, "Not defined for this type of mapping.");
654  return 0;
655  }
656 
658  {
659  ASSERTL0(false, "Not defined for this type of mapping.");
660  static Array<OneD, const int> result;
661  return result;
662  }
663 
664  boost::shared_ptr<AssemblyMap> AssemblyMap::v_LinearSpaceMap(
665  const ExpList &locexp, GlobalSysSolnType solnType)
666  {
667  ASSERTL0(false, "Not defined for this sub class");
668  static boost::shared_ptr<AssemblyMap> result;
669  return result;
670  }
671 
673  {
674  return m_comm;
675  }
676 
677  size_t AssemblyMap::GetHash() const
678  {
679  return m_hash;
680  }
681 
682  int AssemblyMap::GetLocalToGlobalMap(const int i) const
683  {
684  return v_GetLocalToGlobalMap(i);
685  }
686 
688  {
689  return v_GetGlobalToUniversalMap(i);
690  }
691 
693  {
695  }
696 
698  {
699  return v_GetLocalToGlobalMap();
700  }
701 
703  {
704  return v_GetGlobalToUniversalMap();
705  }
706 
708  {
710  }
711 
713  {
714  return v_GetLocalToGlobalSign(i);
715  }
716 
718  {
719  return v_GetLocalToGlobalSign();
720  }
721 
723  const Array<OneD, const NekDouble>& loc,
724  Array<OneD, NekDouble>& global) const
725  {
726  v_LocalToGlobal(loc,global);
727  }
728 
730  const NekVector<NekDouble>& loc,
731  NekVector< NekDouble>& global) const
732  {
733  v_LocalToGlobal(loc,global);
734  }
735 
737  const Array<OneD, const NekDouble>& global,
738  Array<OneD, NekDouble>& loc) const
739  {
740  v_GlobalToLocal(global,loc);
741  }
742 
744  const NekVector<NekDouble>& global,
745  NekVector< NekDouble>& loc) const
746  {
747  v_GlobalToLocal(global,loc);
748  }
749 
751  const Array<OneD, const NekDouble> &loc,
752  Array<OneD, NekDouble> &global) const
753  {
754  v_Assemble(loc,global);
755  }
756 
758  const NekVector<NekDouble>& loc,
759  NekVector< NekDouble>& global) const
760  {
761  v_Assemble(loc,global);
762  }
763 
765  Array<OneD, NekDouble>& pGlobal) const
766  {
767  v_UniversalAssemble(pGlobal);
768  }
769 
771  NekVector< NekDouble>& pGlobal) const
772  {
773  v_UniversalAssemble(pGlobal);
774  }
775 
777  Array<OneD, NekDouble>& pGlobal,
778  int offset) const
779  {
780  v_UniversalAssemble(pGlobal, offset);
781  }
782 
784  {
785  return v_GetFullSystemBandWidth();
786  }
787 
789  {
790  return v_GetNumNonDirVertexModes();
791  }
792 
794  {
795  return v_GetNumNonDirEdgeModes();
796  }
797 
799  {
800  return v_GetNumNonDirFaceModes();
801  }
802 
804  {
805  return v_GetNumDirEdges();
806  }
807 
809  {
810  return v_GetNumDirFaces();
811  }
812 
814  {
815  return v_GetNumNonDirEdges();
816  }
817 
819  {
820  return v_GetNumNonDirFaces();
821  }
822 
824  {
825  return v_GetExtraDirEdges();
826  }
827 
828  boost::shared_ptr<AssemblyMap> AssemblyMap::LinearSpaceMap(const ExpList &locexp, GlobalSysSolnType solnType)
829  {
830  return v_LinearSpaceMap(locexp, solnType);
831  }
832 
833  int AssemblyMap::GetLocalToGlobalBndMap(const int i) const
834  {
835  return m_localToGlobalBndMap[i];
836  }
837 
838  const Array<OneD,const int>&
840  {
841  return m_localToGlobalBndMap;
842  }
843 
845  {
846  return m_signChange;
847  }
848 
849 
852  {
853  return m_localToGlobalBndSign;
854  }
855 
857  {
859  }
860 
862  {
864  }
865 
867  {
868  if(m_signChange)
869  {
870  return m_localToGlobalBndSign[i];
871  }
872  else
873  {
874  return 1.0;
875  }
876  }
877 
878 
880  const int i)
881  {
883  }
884 
885 
887  const int i)
888  {
889  ASSERTL1(i < m_bndCondTraceToGlobalTraceMap.num_elements(),
890  "Index out of range.");
892  }
893 
896  {
897  return m_bndCondTraceToGlobalTraceMap;
898  }
899 
901  {
902  if(m_signChange)
903  {
905  }
906  else
907  {
908  return 1.0;
909  }
910  }
911 
912  const Array<OneD,const int>&
914  {
916  }
917 
918 
920  {
922  }
923 
924 
926  {
927  return m_numLocalDirBndCoeffs;
928  }
929 
931  {
932  return m_numLocalBndCoeffs;
933  }
934 
936  {
937  return m_numGlobalBndCoeffs;
938  }
939 
941  {
942  return m_numLocalCoeffs;
943  }
944 
946  {
947  return m_numGlobalCoeffs;
948  }
949 
951  {
952  return m_systemSingular;
953  }
954 
956  const NekVector<NekDouble>& global,
958  int offset) const
959  {
960  GlobalToLocalBnd(global.GetPtr(), loc.GetPtr(), offset);
961  }
962 
963 
965  const NekVector<NekDouble>& global,
966  NekVector<NekDouble>& loc) const
967  {
968  GlobalToLocalBnd(global.GetPtr(), loc.GetPtr());
969  }
970 
971 
973  const Array<OneD, const NekDouble>& global,
974  Array<OneD,NekDouble>& loc, int offset) const
975  {
976  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
977  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs-offset,"Global vector is not of correct dimension");
978 
979  // offset input data by length "offset" for Dirichlet boundary conditions.
981  Vmath::Vcopy(m_numGlobalBndCoeffs-offset, global.get(), 1, tmp.get() + offset, 1);
982 
983  if(m_signChange)
984  {
986  }
987  else
988  {
989  Vmath::Gathr(m_numLocalBndCoeffs, tmp.get(), m_localToGlobalBndMap.get(), loc.get());
990  }
991  }
992 
993 
995  const Array<OneD, const NekDouble>& global,
996  Array<OneD,NekDouble>& loc) const
997  {
998  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
999  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs,"Global vector is not of correct dimension");
1000 
1001  if(m_signChange)
1002  {
1003  Vmath::Gathr(m_numLocalBndCoeffs, m_localToGlobalBndSign.get(), global.get(), m_localToGlobalBndMap.get(), loc.get());
1004  }
1005  else
1006  {
1007  Vmath::Gathr(m_numLocalBndCoeffs, global.get(), m_localToGlobalBndMap.get(), loc.get());
1008  }
1009  }
1010 
1012  const NekVector<NekDouble>& loc,
1013  NekVector<NekDouble>& global,
1014  int offset) const
1015  {
1016  LocalBndToGlobal(loc.GetPtr(), global.GetPtr(), offset);
1017  }
1018 
1019 
1021  const Array<OneD, const NekDouble>& loc,
1022  Array<OneD,NekDouble>& global,
1023  int offset) const
1024  {
1025  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1026  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs-offset,"Global vector is not of correct dimension");
1027 
1028  // offset input data by length "offset" for Dirichlet boundary conditions.
1030 
1031  if(m_signChange)
1032  {
1034  }
1035  else
1036  {
1037  Vmath::Scatr(m_numLocalBndCoeffs, loc.get(), m_localToGlobalBndMap.get(), tmp.get());
1038  }
1039 
1040  UniversalAssembleBnd(tmp);
1041  Vmath::Vcopy(m_numGlobalBndCoeffs-offset, tmp.get()+offset, 1, global.get(), 1);
1042  }
1043 
1045  const NekVector<NekDouble>& loc,
1046  NekVector<NekDouble>& global) const
1047  {
1048  LocalBndToGlobal(loc.GetPtr(), global.GetPtr());
1049  }
1050 
1051 
1053  const Array<OneD, const NekDouble>& loc,
1054  Array<OneD,NekDouble>& global) const
1055  {
1056  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1057  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs,"Global vector is not of correct dimension");
1058 
1059  if(m_signChange)
1060  {
1061  Vmath::Scatr(m_numLocalBndCoeffs, m_localToGlobalBndSign.get(), loc.get(), m_localToGlobalBndMap.get(), global.get());
1062  }
1063  else
1064  {
1065  Vmath::Scatr(m_numLocalBndCoeffs, loc.get(), m_localToGlobalBndMap.get(), global.get());
1066  }
1067  }
1068 
1069 
1071  const NekVector<NekDouble>& loc,
1072  NekVector<NekDouble>& global, int offset) const
1073  {
1074  AssembleBnd(loc.GetPtr(), global.GetPtr(), offset);
1075  }
1076 
1077 
1079  const NekVector<NekDouble>& loc,
1080  NekVector<NekDouble>& global) const
1081  {
1082  AssembleBnd(loc.GetPtr(), global.GetPtr());
1083  }
1084 
1085 
1087  const Array<OneD,const NekDouble>& loc,
1088  Array<OneD, NekDouble>& global, int offset) const
1089  {
1090  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local array is not of correct dimension");
1091  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs-offset,"Global array is not of correct dimension");
1093 
1094  if(m_signChange)
1095  {
1097  }
1098  else
1099  {
1100  Vmath::Assmb(m_numLocalBndCoeffs,loc.get(), m_localToGlobalBndMap.get(), tmp.get());
1101  }
1102  UniversalAssembleBnd(tmp);
1103  Vmath::Vcopy(m_numGlobalBndCoeffs-offset, tmp.get() + offset, 1, global.get(), 1);
1104  }
1105 
1106 
1108  const Array<OneD, const NekDouble>& loc,
1109  Array<OneD, NekDouble>& global) const
1110  {
1111  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1112  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs,"Global vector is not of correct dimension");
1113 
1114  Vmath::Zero(m_numGlobalBndCoeffs, global.get(), 1);
1115 
1116  if(m_signChange)
1117  {
1119  loc.get(), m_localToGlobalBndMap.get(), global.get());
1120  }
1121  else
1122  {
1123  Vmath::Assmb(m_numLocalBndCoeffs,loc.get(), m_localToGlobalBndMap.get(), global.get());
1124  }
1125  UniversalAssembleBnd(global);
1126  }
1127 
1129  Array<OneD, NekDouble>& pGlobal) const
1130  {
1131  ASSERTL1(pGlobal.num_elements() == m_numGlobalBndCoeffs,
1132  "Wrong size.");
1133  Gs::Gather(pGlobal, Gs::gs_add, m_bndGsh);
1134  }
1135 
1137  NekVector< NekDouble>& pGlobal) const
1138  {
1139  UniversalAssembleBnd(pGlobal.GetPtr());
1140  }
1141 
1143  Array<OneD, NekDouble>& pGlobal,
1144  int offset) const
1145  {
1146  Array<OneD, NekDouble> tmp(offset);
1147  if (offset > 0) Vmath::Vcopy(offset, pGlobal, 1, tmp, 1);
1148  UniversalAssembleBnd(pGlobal);
1149  if (offset > 0) Vmath::Vcopy(offset, tmp, 1, pGlobal, 1);
1150  }
1151 
1153  {
1154  return m_bndSystemBandWidth;
1155  }
1156 
1158  {
1159  return m_staticCondLevel;
1160  }
1161 
1163  {
1164  return m_numPatches;
1165  }
1166 
1169  {
1171  }
1172 
1173 
1176  {
1178  }
1179 
1180  const AssemblyMapSharedPtr
1182  {
1184  }
1185 
1186  const PatchMapSharedPtr&
1188  const
1189  {
1190  return m_patchMapFromPrevLevel;
1191  }
1192 
1194  {
1195  return !( m_nextLevelLocalToGlobalMap.get() );
1196  }
1197 
1198 
1200  {
1201  return m_solnType;
1202  }
1203 
1205  {
1206  return m_preconType;
1207  }
1208 
1210  {
1211  return m_iterativeTolerance;
1212  }
1213 
1215  {
1216  return m_maxIterations;
1217  }
1218 
1220  {
1221  return m_successiveRHS;
1222  }
1223 
1225  const Array<OneD, const NekDouble>& global,
1226  Array<OneD,NekDouble>& loc)
1227  {
1228  ASSERTL1(loc.num_elements() >= m_numLocalBndCoeffs,"Local vector is not of correct dimension");
1229  ASSERTL1(global.num_elements() >= m_numGlobalBndCoeffs,"Global vector is not of correct dimension");
1230 
1231  Vmath::Gathr(m_numLocalBndCoeffs, global.get(), m_localToGlobalBndMap.get(), loc.get());
1232  }
1233 
1235  std::ostream &out, std::string variable) const
1236  {
1238  = m_session->GetComm()->GetRowComm();
1239  bool isRoot = vRowComm->GetRank() == 0;
1240  int n = vRowComm->GetSize();
1241  int i;
1242 
1243  // Determine number of global degrees of freedom.
1244  int globBndCnt = 0, globDirCnt = 0;
1245 
1246  for (i = 0; i < m_numGlobalBndCoeffs; ++i)
1247  {
1249  {
1250  globBndCnt++;
1251 
1252  if (i < m_numGlobalDirBndCoeffs)
1253  {
1254  globDirCnt++;
1255  }
1256  }
1257  }
1258 
1259  int globCnt = m_numGlobalCoeffs - m_numGlobalBndCoeffs + globBndCnt;
1260 
1261  // Calculate maximum valency
1263  Array<OneD, NekDouble> tmpGlob(m_numGlobalBndCoeffs, 0.0);
1264 
1265  Vmath::Assmb(m_numLocalBndCoeffs, tmpLoc.get(), m_localToGlobalBndMap.get(), tmpGlob.get());
1266  UniversalAssembleBnd(tmpGlob);
1267 
1268  int totGlobDof = globCnt;
1269  int totGlobBndDof = globBndCnt;
1270  int totGlobDirDof = globDirCnt;
1271  int totLocalDof = m_numLocalCoeffs;
1272  int totLocalBndDof = m_numLocalBndCoeffs;
1273  int totLocalDirDof = m_numLocalDirBndCoeffs;
1274 
1275  int meanValence = 0;
1276  int maxValence = 0;
1277  int minValence = 10000000;
1278  for (int i = 0; i < m_numGlobalBndCoeffs; ++i)
1279  {
1281  {
1282  continue;
1283  }
1284 
1285  if (tmpGlob[i] > maxValence)
1286  {
1287  maxValence = tmpGlob[i];
1288  }
1289  if (tmpGlob[i] < minValence)
1290  {
1291  minValence = tmpGlob[i];
1292  }
1293  meanValence += tmpGlob[i];
1294  }
1295 
1296  vRowComm->AllReduce(maxValence, LibUtilities::ReduceMax);
1297  vRowComm->AllReduce(minValence, LibUtilities::ReduceMin);
1298  vRowComm->AllReduce(meanValence, LibUtilities::ReduceSum);
1299  vRowComm->AllReduce(totGlobDof, LibUtilities::ReduceSum);
1300  vRowComm->AllReduce(totGlobBndDof, LibUtilities::ReduceSum);
1301  vRowComm->AllReduce(totGlobDirDof, LibUtilities::ReduceSum);
1302  vRowComm->AllReduce(totLocalDof, LibUtilities::ReduceSum);
1303  vRowComm->AllReduce(totLocalBndDof, LibUtilities::ReduceSum);
1304  vRowComm->AllReduce(totLocalDirDof, LibUtilities::ReduceSum);
1305 
1306  meanValence /= totGlobBndDof;
1307 
1308  if (isRoot)
1309  {
1310  out << "Assembly map statistics for field " << variable << ":"
1311  << endl;
1312  out << " - Number of local/global dof : "
1313  << totLocalDof << " " << totGlobDof << endl;
1314  out << " - Number of local/global boundary dof : "
1315  << totLocalBndDof << " " << totGlobBndDof << endl;
1316  out << " - Number of local/global Dirichlet dof : "
1317  << totLocalDirDof << " " << totGlobDirDof << endl;
1318  out << " - dof valency (min/max/mean) : "
1319  << minValence << " " << maxValence << " " << meanValence
1320  << endl;
1321 
1322  if (n > 1)
1323  {
1324  NekDouble mean = m_numLocalCoeffs, mean2 = mean * mean;
1325  NekDouble minval = mean, maxval = mean;
1326  Array<OneD, NekDouble> tmp(1);
1327 
1328  for (i = 1; i < n; ++i)
1329  {
1330  vRowComm->Recv(i, tmp);
1331  mean += tmp[0];
1332  mean2 += tmp[0]*tmp[0];
1333 
1334  if (tmp[0] > maxval)
1335  {
1336  maxval = tmp[0];
1337  }
1338  if (tmp[0] < minval)
1339  {
1340  minval = tmp[0];
1341  }
1342  }
1343 
1344  out << " - Local dof dist. (min/max/mean/dev) : "
1345  << minval << " " << maxval << " " << (mean / n) << " "
1346  << sqrt(mean2/n - mean*mean/n/n) << endl;
1347 
1348  vRowComm->Block();
1349 
1350  mean = minval = maxval = m_numLocalBndCoeffs;
1351  mean2 = mean * mean;
1352 
1353  for (i = 1; i < n; ++i)
1354  {
1355  vRowComm->Recv(i, tmp);
1356  mean += tmp[0];
1357  mean2 += tmp[0]*tmp[0];
1358 
1359  if (tmp[0] > maxval)
1360  {
1361  maxval = tmp[0];
1362  }
1363  if (tmp[0] < minval)
1364  {
1365  minval = tmp[0];
1366  }
1367  }
1368 
1369  out << " - Local bnd dof dist. (min/max/mean/dev) : "
1370  << minval << " " << maxval << " " << (mean / n) << " "
1371  << sqrt(mean2/n - mean*mean/n/n) << endl;
1372  }
1373  }
1374  else
1375  {
1376  Array<OneD, NekDouble> tmp(1);
1377  tmp[0] = m_numLocalCoeffs;
1378  vRowComm->Send(0, tmp);
1379  vRowComm->Block();
1380  tmp[0] = m_numLocalBndCoeffs;
1381  vRowComm->Send(0, tmp);
1382  }
1383  }
1384  } // namespace
1385 } // namespace
const Array< OneD, const int > & GetBndCondTraceToGlobalTraceMap()
PreconditionerType GetPreconType() const
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
bool m_systemSingular
Flag indicating if the system is singular or not.
Definition: AssemblyMap.h:320
const Array< OneD, const int > & GetGlobalToUniversalMap()
bool m_signChange
Flag indicating if modes require sign reversal.
Definition: AssemblyMap.h:345
void PrintStats(std::ostream &out, std::string variable) const
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:368
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:314
LibUtilities::CommSharedPtr GetComm()
Retrieves the communicator.
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: AssemblyMap.h:306
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:224
void Gathr(int n, const T *x, const int *y, T *z)
Gather vector z[i] = x[y[i]].
Definition: Vmath.cpp:630
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:22
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:356
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:331
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:408
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:395
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:377
size_t m_hash
Hash for map.
Definition: AssemblyMap.h:309
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:53
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:365
NekDouble m_iterativeTolerance
Tolerance for iterative solver.
Definition: AssemblyMap.h:374
Array< OneD, unsigned int > m_numLocalBndCoeffsPerPatch
The number of bnd dofs per patch.
Definition: AssemblyMap.h:388
boost::shared_ptr< AssemblyMap > LinearSpaceMap(const ExpList &locexp, GlobalSysSolnType solnType)
GlobalSysSolnType m_solnType
The solution type of the global system.
Definition: AssemblyMap.h:363
Array< OneD, NekDouble > m_bndCondCoeffsToGlobalCoeffsSign
Integer map of bnd cond coeffs to global coefficients.
Definition: AssemblyMap.h:354
int m_numGlobalDirBndCoeffs
Number of Global Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:318
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:659
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:695
void LocalToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
virtual void v_LocalToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
The number of int dofs per patch.
Definition: AssemblyMap.h:390
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:348
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:316
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:352
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:312
int m_staticCondLevel
The level of recursion in the case of multi-level static condensation.
Definition: AssemblyMap.h:384
Array< OneD, NekDouble > m_localToGlobalBndSign
Integer sign of local boundary coeffs to global space.
Definition: AssemblyMap.h:350
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:371
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:303
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:358
Array< OneD, int > m_globalToUniversalBndMapUnique
Integer map of unique process coeffs to universal space (signed)
Definition: AssemblyMap.h:360
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:342
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:359
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:218
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
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:386
bool GetSingularSystem() const
Retrieves if the system is singular (true) or not (false)