Nektar++
AssemblyMap.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: AssemblyMap.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Assembly (e.g. local to global) base mapping routines
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include <boost/algorithm/string.hpp>
36
38
39using namespace std;
40
42{
43/**
44 * @class AssemblyMap
45 * This class acts as a base class for constructing mappings between
46 * local, global and boundary degrees of freedom. It holds the storage
47 * for the maps and provides the accessors needed to retrieve them.
48 *
49 * There are two derived classes: AssemblyMapCG and
50 * AssemblyMapDG. These perform the actual construction of the
51 * maps within their specific contexts.
52 *
53 */
54
55/// Rounds a double precision number to an integer.
57{
58 return int(x > 0.0 ? x + 0.5 : x - 0.5);
59}
60
61/// Rounds an array of double precision numbers to integers.
63 Array<OneD, int> outarray)
64{
65 int size = inarray.size();
66 ASSERTL1(outarray.size() >= 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 */
80 : m_session(), m_comm(), m_hash(0), m_numLocalBndCoeffs(0),
81 m_numGlobalBndCoeffs(0), m_numLocalDirBndCoeffs(0),
82 m_numGlobalDirBndCoeffs(0), m_solnType(eNoSolnType),
83 m_bndSystemBandWidth(0), m_successiveRHS(0),
84 m_linSysIterSolver("ConjugateGradient"), m_gsh(nullptr), m_bndGsh(nullptr)
85{
86}
87
90 const std::string variable)
91 : m_session(pSession), m_comm(comm), m_variable(variable), m_hash(0),
92 m_numLocalBndCoeffs(0), m_numGlobalBndCoeffs(0),
93 m_numLocalDirBndCoeffs(0), m_numGlobalDirBndCoeffs(0),
94 m_bndSystemBandWidth(0), m_successiveRHS(0),
95 m_linSysIterSolver("ConjugateGradient"), m_gsh(nullptr), m_bndGsh(nullptr)
96{
97 // Default value from Solver Info
99 pSession->GetSolverInfoAsEnum<GlobalSysSolnType>("GlobalSysSoln");
100
101 // Override values with data from GlobalSysSolnInfo section
102 if (pSession->DefinesGlobalSysSolnInfo(variable, "GlobalSysSoln"))
103 {
104 std::string sysSoln =
105 pSession->GetGlobalSysSolnInfo(variable, "GlobalSysSoln");
106 m_solnType = pSession->GetValueAsEnum<GlobalSysSolnType>(
107 "GlobalSysSoln", sysSoln);
108 }
109
110 // Set up preconditioner with SysSoln as override then solverinfo otherwise
111 // set a default as diagonal
112 if (pSession->DefinesGlobalSysSolnInfo(variable, "Preconditioner"))
113 {
115 pSession->GetGlobalSysSolnInfo(variable, "Preconditioner");
116 }
117 else if (pSession->DefinesSolverInfo("Preconditioner"))
118 {
119 m_preconType = pSession->GetSolverInfo("Preconditioner");
120 }
121 else
122 { // Possibly preconditioner is default registered in
123 // GlobLinSysIterative.cpp
124 m_preconType = "Diagonal";
125 }
126
127 if (pSession->DefinesGlobalSysSolnInfo(variable,
128 "IterativeSolverTolerance"))
129 {
130 m_iterativeTolerance = boost::lexical_cast<NekDouble>(
131 pSession->GetGlobalSysSolnInfo(variable, "IterativeSolverTolerance")
132 .c_str());
133 }
134 else
135 {
136 pSession->LoadParameter("IterativeSolverTolerance",
139 }
140
141 if (pSession->DefinesGlobalSysSolnInfo(variable, "AbsoluteTolerance"))
142 {
143 std::string abstol =
144 pSession->GetGlobalSysSolnInfo(variable, "AbsoluteTolerance");
146 boost::iequals(boost::to_upper_copy(abstol), "TRUE");
147 }
148 else
149 {
150 m_isAbsoluteTolerance = false;
151 }
152
153 if (pSession->DefinesGlobalSysSolnInfo(variable, "SuccessiveRHS"))
154 {
155 m_successiveRHS = boost::lexical_cast<int>(
156 pSession->GetGlobalSysSolnInfo(variable, "SuccessiveRHS").c_str());
157 }
158 else
159 {
160 pSession->LoadParameter("SuccessiveRHS", m_successiveRHS, 0);
161 }
162
163 if (pSession->DefinesGlobalSysSolnInfo(variable, "LinSysIterSolver"))
164 {
166 pSession->GetGlobalSysSolnInfo(variable, "LinSysIterSolver");
167 }
168 else if (pSession->DefinesSolverInfo("LinSysIterSolver"))
169 {
170 m_linSysIterSolver = pSession->GetSolverInfo("LinSysIterSolver");
171 }
172 else
173 {
174 m_linSysIterSolver = "ConjugateGradient";
175 }
176}
177
178/**
179 * Create a new level of mapping using the information in
180 * multiLevelGraph and performing the following steps:
181 */
183 AssemblyMap *oldLevelMap,
184 const BottomUpSubStructuredGraphSharedPtr &multiLevelGraph)
185 : m_session(oldLevelMap->m_session), m_comm(oldLevelMap->GetComm()),
186 m_hash(0), m_solnType(oldLevelMap->m_solnType),
187 m_preconType(oldLevelMap->m_preconType),
188 m_iterativeTolerance(oldLevelMap->m_iterativeTolerance),
189 m_successiveRHS(oldLevelMap->m_successiveRHS),
190 m_linSysIterSolver(oldLevelMap->m_linSysIterSolver),
191 m_gsh(oldLevelMap->m_gsh), m_bndGsh(oldLevelMap->m_bndGsh),
192 m_lowestStaticCondLevel(oldLevelMap->m_lowestStaticCondLevel)
193{
194 int i;
195 int j;
196 int cnt;
197
198 //--------------------------------------------------------------
199 // -- Extract information from the input argument
200 int numGlobalBndCoeffsOld = oldLevelMap->GetNumGlobalBndCoeffs();
201 int numGlobalDirBndCoeffsOld = oldLevelMap->GetNumGlobalDirBndCoeffs();
202 int numLocalBndCoeffsOld = oldLevelMap->GetNumLocalBndCoeffs();
203 int numLocalDirBndCoeffsOld = oldLevelMap->GetNumLocalDirBndCoeffs();
204 bool signChangeOld = oldLevelMap->GetSignChange();
205
206 int staticCondLevelOld = oldLevelMap->GetStaticCondLevel();
207 int numPatchesOld = oldLevelMap->GetNumPatches();
208 GlobalSysSolnType solnTypeOld = oldLevelMap->GetGlobalSysSolnType();
209 const Array<OneD, const unsigned int> &numLocalBndCoeffsPerPatchOld =
210 oldLevelMap->GetNumLocalBndCoeffsPerPatch();
211 //--------------------------------------------------------------
212
213 //--------------------------------------------------------------
214 int newLevel = staticCondLevelOld + 1;
215 /** - STEP 1: setup a mask array to determine to which patch
216 * of the new level every patch of the current
217 * level belongs. To do so we make four arrays,
218 * #gloPatchMask, #globHomPatchMask,
219 * #locPatchMask_NekDouble and #locPatchMask.
220 * These arrays are then used to check which local
221 * dofs of the old level belong to which patch of
222 * the new level
223 */
224 Array<OneD, NekDouble> globPatchMask(numGlobalBndCoeffsOld, -1.0);
225 Array<OneD, NekDouble> globHomPatchMask(globPatchMask +
226 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,
244 locPatchMask_NekDouble);
245
246 // Convert the result to an array of integers rather than doubles
247 RoundNekDoubleToInt(locPatchMask_NekDouble, locPatchMask);
248
249 /** - STEP 2: We calculate how many local bnd dofs of the
250 * old level belong to the boundaries of each patch at
251 * the new level. We need this information to set up the
252 * mapping between different levels.
253 */
254
255 // Retrieve the number of patches at the next level
256 int numPatchesWithIntNew =
257 multiLevelGraph->GetNpatchesWithInterior(newLevel);
258 int numPatchesNew = numPatchesWithIntNew;
259
260 // Allocate memory to store the number of local dofs associated to
261 // each of elemental boundaries of these patches
262 std::map<int, int> numLocalBndCoeffsPerPatchNew;
263 for (int i = 0; i < numPatchesNew; i++)
264 {
265 numLocalBndCoeffsPerPatchNew[i] = 0;
266 }
267
268 int minval;
269 int maxval;
270 int curPatch;
271 for (i = cnt = 0; i < numPatchesOld; i++)
272 {
273 // For every patch at the current level, the mask array
274 // locPatchMask should be filled with
275 // - the same (positive) number for each entry (which will
276 // correspond to the patch at the next level it belongs to)
277 // - the same (positive) number for each entry, except some
278 // entries that are -1 (the enties correspond to -1, will be
279 // mapped to the local boundary of the next level patch given
280 // by the positive number)
281 // - -1 for all entries. In this case, we will make an
282 // additional patch only consisting of boundaries at the next
283 // level
284 minval =
285 *min_element(&locPatchMask[cnt],
286 &locPatchMask[cnt] + numLocalBndCoeffsPerPatchOld[i]);
287 maxval =
288 *max_element(&locPatchMask[cnt],
289 &locPatchMask[cnt] + numLocalBndCoeffsPerPatchOld[i]);
290 ASSERTL0((minval == maxval) || (minval == -1),
291 "These values should never be the same");
292
293 if (maxval == -1)
294 {
295 curPatch = numPatchesNew;
296 numLocalBndCoeffsPerPatchNew[curPatch] = 0;
297 numPatchesNew++;
298 }
299 else
300 {
301 curPatch = maxval;
302 }
303
304 for (j = 0; j < numLocalBndCoeffsPerPatchOld[i]; j++)
305 {
306 ASSERTL0((locPatchMask[cnt] == maxval) ||
307 (locPatchMask[cnt] == minval),
308 "These values should never be the same");
309 if (locPatchMask[cnt] == -1)
310 {
311 ++numLocalBndCoeffsPerPatchNew[curPatch];
312 }
313 cnt++;
314 }
315 }
316
317 // Count how many local dofs of the old level are mapped
318 // to the local boundary dofs of the new level
321 m_numPatches = numLocalBndCoeffsPerPatchNew.size();
324 multiLevelGraph->GetNintDofsPerPatch(newLevel, m_numLocalIntCoeffsPerPatch);
325
326 for (int i = 0; i < m_numPatches; i++)
327 {
329 (unsigned int)numLocalBndCoeffsPerPatchNew[i];
333 }
334
335 // Also initialise some more data members
336 m_solnType = solnTypeOld;
341 "This method should only be called for in "
342 "case of multi-level static condensation.");
343 m_staticCondLevel = newLevel;
344 m_signChange = signChangeOld;
345 m_numLocalDirBndCoeffs = numLocalDirBndCoeffsOld;
346 m_numGlobalDirBndCoeffs = numGlobalDirBndCoeffsOld;
348 multiLevelGraph->GetInteriorOffset(newLevel) + m_numGlobalDirBndCoeffs;
350 multiLevelGraph->GetNumGlobalDofs(newLevel) + m_numGlobalDirBndCoeffs;
352
356
357 // local to bnd map is just a copy
358 for (int i = 0; i < m_numLocalBndCoeffs; ++i)
359 {
361 }
362
363 // local to int map is just a copy plus offset
364 for (int i = m_numLocalBndCoeffs; i < m_numLocalCoeffs; ++i)
365 {
367 }
368
369 if (m_signChange)
370 {
372 }
373
375 MemoryManager<PatchMap>::AllocateSharedPtr(numLocalBndCoeffsOld);
376
381
382 // Set up an offset array that denotes the offset of the local
383 // boundary degrees of freedom of the next level
384 Array<OneD, int> numLocalBndCoeffsPerPatchOffset(m_numPatches + 1, 0);
385 for (int i = 1; i < m_numPatches + 1; i++)
386 {
387 numLocalBndCoeffsPerPatchOffset[i] +=
388 numLocalBndCoeffsPerPatchOffset[i - 1] +
389 numLocalBndCoeffsPerPatchNew[i - 1];
390 }
391
392 int additionalPatchCnt = numPatchesWithIntNew;
393 int newid;
394 int blockid;
395 bool isBndDof;
397 Array<OneD, int> bndDofPerPatchCnt(m_numPatches, 0);
398
399 for (i = cnt = 0; i < numPatchesOld; i++)
400 {
401 minval =
402 *min_element(&locPatchMask[cnt],
403 &locPatchMask[cnt] + numLocalBndCoeffsPerPatchOld[i]);
404 maxval =
405 *max_element(&locPatchMask[cnt],
406 &locPatchMask[cnt] + numLocalBndCoeffsPerPatchOld[i]);
407 ASSERTL0((minval == maxval) || (minval == -1),
408 "These values should never be the same");
409
410 if (maxval == -1)
411 {
412 curPatch = additionalPatchCnt;
413 additionalPatchCnt++;
414 }
415 else
416 {
417 curPatch = maxval;
418 }
419
420 for (j = 0; j < numLocalBndCoeffsPerPatchOld[i]; j++)
421 {
422 ASSERTL0((locPatchMask[cnt] == maxval) ||
423 (locPatchMask[cnt] == minval),
424 "These values should never be the same");
425
426 sign = oldLevelMap->GetLocalToGlobalBndSign(cnt);
427
428 if (locPatchMask[cnt] == -1)
429 {
430 newid = numLocalBndCoeffsPerPatchOffset[curPatch];
431
432 m_localToGlobalBndMap[newid] =
433 oldLevelMap->GetLocalToGlobalBndMap(cnt);
434
435 if (m_signChange)
436 {
438 }
439
440 blockid = bndDofPerPatchCnt[curPatch];
441 isBndDof = true;
442
443 numLocalBndCoeffsPerPatchOffset[curPatch]++;
444 bndDofPerPatchCnt[curPatch]++;
445 }
446 else
447 {
448 newid = oldLevelMap->GetLocalToGlobalBndMap(cnt) -
450
451 blockid =
452 oldLevelMap->GetLocalToGlobalBndMap(cnt) -
454 multiLevelGraph->GetInteriorOffset(newLevel, curPatch);
455 isBndDof = false;
456 }
457
458 sign = isBndDof ? 1.0 : sign;
459
460 m_patchMapFromPrevLevel->SetPatchMap(cnt, curPatch, blockid,
461 isBndDof, sign);
462 cnt++;
463 }
464 }
465
466 // set up local to local mapping from previous to new level
469
471
472 // Postprocess the computed information - Update the old
473 // level with the mapping to new level
474 // oldLevelMap->SetLocalBndToNextLevelMap(oldLocalBndToNewLevelMap,oldLocalBndToNewLevelSign);
475
476 // - Construct the next level mapping object
477 if (m_staticCondLevel < (multiLevelGraph->GetNlevels() - 1))
478 {
481 multiLevelGraph);
482 }
483}
484
486{
487}
488
489/**
490 * The bandwidth calculated corresponds to what is referred to as
491 * half-bandwidth. If the elements of the matrix are designated as
492 * a_ij, it corresponds to the maximum value of |i-j| for non-zero
493 * a_ij. As a result, the value also corresponds to the number of
494 * sub- or super-diagonals.
495 *
496 * The bandwith can be calculated elementally as it corresponds to the
497 * maximal elemental bandwith (i.e. the maximal difference in global
498 * DOF index for every element).
499 *
500 * We here calculate the bandwith of the global boundary system (as
501 * used for static condensation).
502 */
504{
505 int i, j;
506 int cnt = 0;
507 int locSize;
508 int maxId;
509 int minId;
510 int bwidth = -1;
511 for (i = 0; i < m_numPatches; ++i)
512 {
513 locSize = m_numLocalBndCoeffsPerPatch[i];
514 maxId = -1;
515 minId = m_numLocalBndCoeffs + 1;
516 for (j = 0; j < locSize; j++)
517 {
519 {
520 if (m_localToGlobalBndMap[cnt + j] > maxId)
521 {
522 maxId = m_localToGlobalBndMap[cnt + j];
523 }
524
525 if (m_localToGlobalBndMap[cnt + j] < minId)
526 {
527 minId = m_localToGlobalBndMap[cnt + j];
528 }
529 }
530 }
531 bwidth = (bwidth > (maxId - minId)) ? bwidth : (maxId - minId);
532
533 cnt += locSize;
534 }
535
536 m_bndSystemBandWidth = bwidth;
537}
538
539int AssemblyMap::v_GetLocalToGlobalMap([[maybe_unused]] const int i) const
540{
541 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
542 return 0;
543}
544
545int AssemblyMap::v_GetGlobalToUniversalMap([[maybe_unused]] const int i) const
546{
547 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
548 return 0;
549}
550
552 [[maybe_unused]] const int i) const
553{
554 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
555 return 0;
556}
557
559{
560 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
561 static Array<OneD, const int> result;
562 return result;
563}
564
566{
567 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
568 static Array<OneD, const int> result;
569 return result;
570}
571
573{
574 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
575 static Array<OneD, const int> result;
576 return result;
577}
578
580 [[maybe_unused]] const int i) const
581{
582 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
583 return 0.0;
584}
585
587{
588 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
589 static Array<OneD, NekDouble> result;
590 return result;
591}
592
594 [[maybe_unused]] const Array<OneD, const NekDouble> &loc,
595 [[maybe_unused]] Array<OneD, NekDouble> &global,
596 [[maybe_unused]] bool useComm) const
597{
598 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
599}
600
602 [[maybe_unused]] const NekVector<NekDouble> &loc,
603 [[maybe_unused]] NekVector<NekDouble> &global,
604 [[maybe_unused]] bool useComm) const
605{
606 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
607}
608
610 [[maybe_unused]] const Array<OneD, const NekDouble> &global,
611 [[maybe_unused]] Array<OneD, NekDouble> &loc) const
612{
613 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
614}
615
617 [[maybe_unused]] const NekVector<NekDouble> &global,
618 [[maybe_unused]] NekVector<NekDouble> &loc) const
619{
620 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
621}
622
624 [[maybe_unused]] const Array<OneD, const NekDouble> &loc,
625 [[maybe_unused]] Array<OneD, NekDouble> &global) const
626{
627 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
628}
629
631 [[maybe_unused]] const NekVector<NekDouble> &loc,
632 [[maybe_unused]] NekVector<NekDouble> &global) const
633{
634 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
635}
636
638 [[maybe_unused]] Array<OneD, NekDouble> &pGlobal) const
639{
640 // Do nothing here since multi-level static condensation uses a
641 // AssemblyMap and thus will call this routine in serial.
642}
643
645 [[maybe_unused]] NekVector<NekDouble> &pGlobal) const
646{
647 // Do nothing here since multi-level static condensation uses a
648 // AssemblyMap and thus will call this routine in serial.
649}
650
652 [[maybe_unused]] Array<OneD, NekDouble> &pGlobal,
653 [[maybe_unused]] int offset) const
654{
655 // Do nothing here since multi-level static condensation uses a
656 // AssemblyMap and thus will call this routine in serial.
657}
658
660{
661 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
662 return 0;
663}
664
666{
667 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
668 return 0;
669}
670
672{
673 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
674 return 0;
675}
676
678{
679 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
680 return 0;
681}
682
684{
685 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
686 return 0;
687}
688
690{
691 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
692 return 0;
693}
694
696{
697 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
698 return 0;
699}
700
702{
703 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
704 return 0;
705}
706
708{
709 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
710 static Array<OneD, const int> result;
711 return result;
712}
713
714std::shared_ptr<AssemblyMap> AssemblyMap::v_LinearSpaceMap(
715 [[maybe_unused]] const ExpList &locexp,
716 [[maybe_unused]] GlobalSysSolnType solnType)
717{
718 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
719 static std::shared_ptr<AssemblyMap> result;
720 return result;
721}
722
724{
725 return m_comm;
726}
727
729{
730 return m_variable;
731}
732
734{
735 return m_hash;
736}
737
739{
740 return v_GetLocalToGlobalMap(i);
741}
742
744{
746}
747
749{
751}
752
754{
755 return v_GetLocalToGlobalMap();
756}
757
759{
761}
762
764{
766}
767
769{
770 return v_GetLocalToGlobalSign(i);
771}
772
774{
775 return v_GetLocalToGlobalSign();
776}
777
780 bool useComm) const
781{
782 v_LocalToGlobal(loc, global, useComm);
783}
784
786 NekVector<NekDouble> &global,
787 bool useComm) const
788{
789 v_LocalToGlobal(loc, global, useComm);
790}
791
794{
795 v_GlobalToLocal(global, loc);
796}
797
800{
801 v_GlobalToLocal(global, loc);
802}
803
805 Array<OneD, NekDouble> &global) const
806{
807 v_Assemble(loc, global);
808}
809
811 NekVector<NekDouble> &global) const
812{
813 v_Assemble(loc, global);
814}
815
817{
818 v_UniversalAssemble(pGlobal);
819}
820
822{
823 v_UniversalAssemble(pGlobal);
824}
825
827 int offset) const
828{
829 v_UniversalAssemble(pGlobal, offset);
830}
831
833 Array<OneD, NekDouble> &global) const
834{
836
837 Array<OneD, const int> map = m_patchMapFromPrevLevel->GetNewLevelMap();
839
840 if (global.data() == loc.data())
841 {
842 local = Array<OneD, NekDouble>(map.size(), loc.data());
843 }
844 else
845 {
846 local = loc; // create reference
847 }
848
849 Vmath::Scatr(map.size(), sign.get(), local.get(), map.get(), global.get());
850}
851
854{
856
857 Array<OneD, const int> map = m_patchMapFromPrevLevel->GetNewLevelMap();
859
860 if (global.data() == loc.data())
861 {
862 glo = Array<OneD, NekDouble>(global.size(), global.data());
863 }
864 else
865 {
866 glo = global; // create reference
867 }
868
869 Vmath::Gathr(map.size(), sign.get(), glo.get(), map.get(), loc.get());
870}
871
873 Array<OneD, NekDouble> &global) const
874{
876 Array<OneD, const int> map = m_patchMapFromPrevLevel->GetNewLevelMap();
878
879 if (global.data() == loc.data())
880 {
881 local = Array<OneD, NekDouble>(map.size(), loc.data());
882 }
883 else
884 {
885 local = loc; // create reference
886 }
887
888 // since we are calling mapping from level down from array
889 // the m_numLocaBndCoeffs represents the size of the
890 // boundary elements we need to assemble into
891 Vmath::Zero(m_numLocalCoeffs, global.get(), 1);
892
893 Vmath::Assmb(map.size(), sign.get(), local.get(), map.get(), global.get());
894}
895
897{
899}
900
902{
904}
905
907{
909}
910
912{
914}
915
917{
918 return v_GetNumDirEdges();
919}
920
922{
923 return v_GetNumDirFaces();
924}
925
927{
928 return v_GetNumNonDirEdges();
929}
930
932{
933 return v_GetNumNonDirFaces();
934}
935
937{
938 return v_GetExtraDirEdges();
939}
940
941std::shared_ptr<AssemblyMap> AssemblyMap::LinearSpaceMap(
942 const ExpList &locexp, GlobalSysSolnType solnType)
943{
944 return v_LinearSpaceMap(locexp, solnType);
945}
946
948{
949 return m_localToGlobalBndMap[i];
950}
951
953{
955}
956
958{
959 return m_signChange;
960}
961
963{
965}
966
968{
970}
971
973{
975}
976
978{
979 if (m_signChange)
980 {
981 return m_localToGlobalBndSign[i];
982 }
983 else
984 {
985 return 1.0;
986 }
987}
988
990{
992}
993
995{
997}
998
1000{
1002}
1003
1005{
1006 ASSERTL1(i < m_bndCondIDToGlobalTraceID.size(), "Index out of range.");
1008}
1009
1011{
1013}
1014
1016{
1018}
1019
1021{
1023}
1024
1026{
1027 return m_numLocalBndCoeffs;
1028}
1029
1031{
1032 return m_numGlobalBndCoeffs;
1033}
1034
1036{
1037 return m_numLocalCoeffs;
1038}
1039
1041{
1042 return m_numGlobalCoeffs;
1043}
1044
1046{
1047 return m_systemSingular;
1048}
1049
1051 NekVector<NekDouble> &loc, int offset) const
1052{
1053 GlobalToLocalBnd(global.GetPtr(), loc.GetPtr(), offset);
1054}
1055
1058{
1059 GlobalToLocalBnd(global.GetPtr(), loc.GetPtr());
1060}
1061
1064 int offset) const
1065{
1067 "Local vector is not of correct dimension");
1068 ASSERTL1(global.size() >= m_numGlobalBndCoeffs - offset,
1069 "Global vector is not of correct dimension");
1070
1071 // offset input data by length "offset" for Dirichlet boundary conditions.
1073 Vmath::Vcopy(m_numGlobalBndCoeffs - offset, global.get(), 1,
1074 tmp.get() + offset, 1);
1075
1076 if (m_signChange)
1077 {
1079 tmp.get(), m_localToGlobalBndMap.get(), loc.get());
1080 }
1081 else
1082 {
1084 m_localToGlobalBndMap.get(), loc.get());
1085 }
1086}
1087
1090{
1092 "Local vector is not of correct dimension");
1093 ASSERTL1(global.size() >= m_numGlobalBndCoeffs,
1094 "Global vector is not of correct dimension");
1095
1097 if (global.data() == loc.data())
1098 {
1099 glo = Array<OneD, NekDouble>(m_numLocalBndCoeffs, global.data());
1100 }
1101 else
1102 {
1103 glo = global; // create reference
1104 }
1105
1106 if (m_signChange)
1107 {
1109 glo.get(), m_localToGlobalBndMap.get(), loc.get());
1110 }
1111 else
1112 {
1114 m_localToGlobalBndMap.get(), loc.get());
1115 }
1116}
1117
1119 Array<OneD, NekDouble> &global, int offset,
1120 bool UseComm) const
1121{
1123 "Local vector is not of correct dimension");
1124 ASSERTL1(global.size() >= m_numGlobalBndCoeffs - offset,
1125 "Global vector is not of correct dimension");
1126
1127 // offset input data by length "offset" for Dirichlet boundary conditions.
1129
1130 if (m_signChange)
1131 {
1133 loc.get(), m_localToGlobalBndMap.get(), tmp.get());
1134 }
1135 else
1136 {
1138 m_localToGlobalBndMap.get(), tmp.get());
1139 }
1140
1141 // Ensure each processor has unique value with a max gather.
1142 if (UseComm)
1143 {
1145 }
1146 Vmath::Vcopy(m_numGlobalBndCoeffs - offset, tmp.get() + offset, 1,
1147 global.get(), 1);
1148}
1149
1151 Array<OneD, NekDouble> &global,
1152 bool UseComm) const
1153{
1155 "Local vector is not of correct dimension");
1156 ASSERTL1(global.size() >= m_numGlobalBndCoeffs,
1157 "Global vector is not of correct dimension");
1158
1159 if (m_signChange)
1160 {
1162 loc.get(), m_localToGlobalBndMap.get(), global.get());
1163 }
1164 else
1165 {
1167 m_localToGlobalBndMap.get(), global.get());
1168 }
1169 if (UseComm)
1170 {
1171 Gs::Gather(global, Gs::gs_max, m_bndGsh);
1172 }
1173}
1174
1176 Array<OneD, NekDouble> &locbnd) const
1177{
1178 ASSERTL1(locbnd.size() >= m_numLocalBndCoeffs,
1179 "LocBnd vector is not of correct dimension");
1180 ASSERTL1(local.size() >= m_numLocalCoeffs,
1181 "Local vector is not of correct dimension");
1182
1184 locbnd.get());
1185}
1186
1188 Array<OneD, NekDouble> &locint) const
1189{
1191 "Locint vector is not of correct dimension");
1192 ASSERTL1(local.size() >= m_numLocalCoeffs,
1193 "Local vector is not of correct dimension");
1194
1196 m_localToLocalIntMap.get(), locint.get());
1197}
1198
1200 Array<OneD, NekDouble> &local) const
1201{
1202 ASSERTL1(locbnd.size() >= m_numLocalBndCoeffs,
1203 "LocBnd vector is not of correct dimension");
1204 ASSERTL1(local.size() >= m_numLocalCoeffs,
1205 "Local vector is not of correct dimension");
1206
1208 local.get());
1209}
1210
1212 Array<OneD, NekDouble> &local) const
1213{
1215 "LocBnd vector is not of correct dimension");
1216 ASSERTL1(local.size() >= m_numLocalCoeffs,
1217 "Local vector is not of correct dimension");
1218
1220 m_localToLocalIntMap.get(), local.get());
1221}
1222
1224 NekVector<NekDouble> &global, int offset) const
1225{
1226 AssembleBnd(loc.GetPtr(), global.GetPtr(), offset);
1227}
1228
1230 NekVector<NekDouble> &global) const
1231{
1232 AssembleBnd(loc.GetPtr(), global.GetPtr());
1233}
1234
1236 Array<OneD, NekDouble> &global, int offset) const
1237{
1239 "Local array is not of correct dimension");
1240 ASSERTL1(global.size() >= m_numGlobalBndCoeffs - offset,
1241 "Global array is not of correct dimension");
1243
1244 if (m_signChange)
1245 {
1247 loc.get(), m_localToGlobalBndMap.get(), tmp.get());
1248 }
1249 else
1250 {
1252 m_localToGlobalBndMap.get(), tmp.get());
1253 }
1255 Vmath::Vcopy(m_numGlobalBndCoeffs - offset, tmp.get() + offset, 1,
1256 global.get(), 1);
1257}
1258
1260 Array<OneD, NekDouble> &global) const
1261{
1263 "Local vector is not of correct dimension");
1264 ASSERTL1(global.size() >= m_numGlobalBndCoeffs,
1265 "Global vector is not of correct dimension");
1266
1267 Vmath::Zero(m_numGlobalBndCoeffs, global.get(), 1);
1268
1269 if (m_signChange)
1270 {
1272 loc.get(), m_localToGlobalBndMap.get(), global.get());
1273 }
1274 else
1275 {
1277 m_localToGlobalBndMap.get(), global.get());
1278 }
1279 UniversalAssembleBnd(global);
1280}
1281
1283{
1284 ASSERTL1(pGlobal.size() >= m_numGlobalBndCoeffs, "Wrong size.");
1285 Gs::Gather(pGlobal, Gs::gs_add, m_bndGsh);
1286}
1287
1289{
1290 UniversalAssembleBnd(pGlobal.GetPtr());
1291}
1292
1294 int offset) const
1295{
1296 Array<OneD, NekDouble> tmp(offset);
1297 if (offset > 0)
1298 {
1299 Vmath::Vcopy(offset, pGlobal, 1, tmp, 1);
1300 }
1301 UniversalAssembleBnd(pGlobal);
1302 if (offset > 0)
1303 {
1304 Vmath::Vcopy(offset, tmp, 1, pGlobal, 1);
1305 }
1306}
1307
1309{
1311}
1312
1314{
1315 return m_bndSystemBandWidth;
1316}
1317
1319{
1320 return m_staticCondLevel;
1321}
1322
1324{
1325 return m_numPatches;
1326}
1327
1330{
1332}
1333
1336{
1338}
1339
1341{
1343}
1344
1346{
1348}
1349
1351{
1352 return !(m_nextLevelLocalToGlobalMap.get());
1353}
1354
1356{
1357 return m_solnType;
1358}
1359
1361{
1362 return m_preconType;
1363}
1364
1366{
1367 return m_iterativeTolerance;
1368}
1369
1371{
1372 return m_isAbsoluteTolerance;
1373}
1374
1376{
1377 return m_successiveRHS;
1378}
1379
1381{
1382 return m_linSysIterSolver;
1383}
1384
1387{
1389 "Local vector is not of correct dimension");
1390 ASSERTL1(global.size() >= m_numGlobalBndCoeffs,
1391 "Global vector is not of correct dimension");
1392
1394 loc.get());
1395}
1396
1397void AssemblyMap::PrintStats(std::ostream &out, std::string variable,
1398 bool printHeader) const
1399{
1400 LibUtilities::CommSharedPtr vRowComm = m_session->GetComm()->GetRowComm();
1401 bool isRoot = m_session->GetComm()->IsParallelInTime()
1402 ? m_session->GetComm()->GetRank() == 0
1403 : vRowComm->GetRank() == 0;
1404 int n = vRowComm->GetSize();
1405 int i;
1406
1407 // Determine number of global degrees of freedom.
1408 int globBndCnt = 0, globDirCnt = 0;
1409
1410 for (i = 0; i < m_numGlobalBndCoeffs; ++i)
1411 {
1413 {
1414 globBndCnt++;
1415
1417 {
1418 globDirCnt++;
1419 }
1420 }
1421 }
1422
1423 int globCnt = m_numGlobalCoeffs - m_numGlobalBndCoeffs + globBndCnt;
1424
1425 // Calculate maximum valency
1428
1430 tmpGlob.get());
1431 UniversalAssembleBnd(tmpGlob);
1432
1433 int totGlobDof = globCnt;
1434 int totGlobBndDof = globBndCnt;
1435 int totGlobDirDof = globDirCnt;
1436 int totLocalDof = m_numLocalCoeffs;
1437 int totLocalBndDof = m_numLocalBndCoeffs;
1438 int totLocalDirDof = m_numLocalDirBndCoeffs;
1439
1440 int meanValence = 0;
1441 int maxValence = 0;
1442 int minValence = 10000000;
1443 for (int i = 0; i < m_numGlobalBndCoeffs; ++i)
1444 {
1446 {
1447 continue;
1448 }
1449
1450 if (tmpGlob[i] > maxValence)
1451 {
1452 maxValence = tmpGlob[i];
1453 }
1454 if (tmpGlob[i] < minValence)
1455 {
1456 minValence = tmpGlob[i];
1457 }
1458 meanValence += tmpGlob[i];
1459 }
1460
1461 vRowComm->AllReduce(maxValence, LibUtilities::ReduceMax);
1462 vRowComm->AllReduce(minValence, LibUtilities::ReduceMin);
1463 vRowComm->AllReduce(meanValence, LibUtilities::ReduceSum);
1464 vRowComm->AllReduce(totGlobDof, LibUtilities::ReduceSum);
1465 vRowComm->AllReduce(totGlobBndDof, LibUtilities::ReduceSum);
1466 vRowComm->AllReduce(totGlobDirDof, LibUtilities::ReduceSum);
1467 vRowComm->AllReduce(totLocalDof, LibUtilities::ReduceSum);
1468 vRowComm->AllReduce(totLocalBndDof, LibUtilities::ReduceSum);
1469 vRowComm->AllReduce(totLocalDirDof, LibUtilities::ReduceSum);
1470
1471 meanValence /= totGlobBndDof;
1472
1473 if (isRoot)
1474 {
1475 if (printHeader)
1476 {
1477 out << "Assembly map statistics for field " << variable << ":"
1478 << endl;
1479 }
1480
1481 out << " - Number of local/global dof : " << totLocalDof
1482 << " " << totGlobDof << endl;
1483 out << " - Number of local/global boundary dof : " << totLocalBndDof
1484 << " " << totGlobBndDof << endl;
1485 out << " - Number of local/global Dirichlet dof : " << totLocalDirDof
1486 << " " << totGlobDirDof << endl;
1487 out << " - dof valency (min/max/mean) : " << minValence
1488 << " " << maxValence << " " << meanValence << endl;
1489
1490 if (n > 1)
1491 {
1492 NekDouble mean = m_numLocalCoeffs, mean2 = mean * mean;
1493 NekDouble minval = mean, maxval = mean;
1495
1496 for (i = 1; i < n; ++i)
1497 {
1498 vRowComm->Recv(i, tmp);
1499 mean += tmp[0];
1500 mean2 += tmp[0] * tmp[0];
1501
1502 if (tmp[0] > maxval)
1503 {
1504 maxval = tmp[0];
1505 }
1506 if (tmp[0] < minval)
1507 {
1508 minval = tmp[0];
1509 }
1510 }
1511
1512 if (maxval > 0.1)
1513 {
1514 out << " - Local dof dist. (min/max/mean/dev) : " << minval
1515 << " " << maxval << " " << (mean / n) << " "
1516 << sqrt(mean2 / n - mean * mean / n / n) << endl;
1517 }
1518
1519 vRowComm->Block();
1520
1521 mean = minval = maxval = m_numLocalBndCoeffs;
1522 mean2 = mean * mean;
1523
1524 for (i = 1; i < n; ++i)
1525 {
1526 vRowComm->Recv(i, tmp);
1527 mean += tmp[0];
1528 mean2 += tmp[0] * tmp[0];
1529
1530 if (tmp[0] > maxval)
1531 {
1532 maxval = tmp[0];
1533 }
1534 if (tmp[0] < minval)
1535 {
1536 minval = tmp[0];
1537 }
1538 }
1539
1540 out << " - Local bnd dof dist. (min/max/mean/dev) : " << minval
1541 << " " << maxval << " " << (mean / n) << " "
1542 << sqrt(mean2 / n - mean * mean / n / n) << endl;
1543 }
1544 }
1545 else
1546 {
1548 tmp[0] = m_numLocalCoeffs;
1549 vRowComm->Send(0, tmp);
1550 vRowComm->Block();
1551 tmp[0] = m_numLocalBndCoeffs;
1552 vRowComm->Send(0, tmp);
1553 }
1554
1555 // Either we have no more levels in the static condensation, or we
1556 // are not multi-level.
1558 {
1559 return;
1560 }
1561
1562 int level = 2;
1564 while (tmp->m_nextLevelLocalToGlobalMap)
1565 {
1566 tmp = tmp->m_nextLevelLocalToGlobalMap;
1567 ++level;
1568 }
1569
1570 // Print out multi-level static condensation information.
1571 if (n > 1)
1572 {
1573 if (isRoot)
1574 {
1575 NekDouble mean = level, mean2 = mean * mean;
1576 int minval = level, maxval = level;
1577
1578 Array<OneD, NekDouble> tmpRecv(1);
1579 for (i = 1; i < n; ++i)
1580 {
1581 vRowComm->Recv(i, tmpRecv);
1582 mean += tmpRecv[0];
1583 mean2 += tmpRecv[0] * tmpRecv[0];
1584
1585 if (tmpRecv[0] > maxval)
1586 {
1587 maxval = (int)(tmpRecv[0] + 0.5);
1588 }
1589 if (tmpRecv[0] < minval)
1590 {
1591 minval = (int)(tmpRecv[0] + 0.5);
1592 }
1593 }
1594
1595 out << " - M-level sc. dist. (min/max/mean/dev) : " << minval
1596 << " " << maxval << " " << (mean / n) << " "
1597 << sqrt(mean2 / n - mean * mean / n / n) << endl;
1598 }
1599 else
1600 {
1601 Array<OneD, NekDouble> tmpSend(1);
1602 tmpSend[0] = level;
1603 vRowComm->Send(0, tmpSend);
1604 }
1605 }
1606 else
1607 {
1608 if (isRoot)
1609 {
1610 out << " - Number of static cond. levels : " << level
1611 << endl;
1612 }
1613 }
1614
1615 if (isRoot)
1616 {
1617 out << "Stats at lowest static cond. level:" << endl;
1618 }
1619 tmp->PrintStats(out, variable, false);
1620}
1621} // namespace Nektar::MultiRegions
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:47
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
Base class for constructing local to global mapping of degrees of freedom.
Definition: AssemblyMap.h:56
void PatchLocalToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
const Array< OneD, const int > & GetExtraDirEdges()
GlobalSysSolnType m_solnType
The solution type of the global system.
Definition: AssemblyMap.h:402
int GetNumPatches() const
Returns the number of patches in this static condensation level.
const Array< OneD, const int > & GetGlobalToUniversalMapUnique()
int GetNumGlobalCoeffs() const
Returns the total number of global coefficients.
int m_numLocalCoeffs
Total number of local coefficients.
Definition: AssemblyMap.h:364
Array< OneD, int > m_bndCondCoeffsToLocalTraceMap
Integer map of bnd cond coeff to local trace coeff.
Definition: AssemblyMap.h:393
const Array< OneD, const unsigned int > & GetNumLocalBndCoeffsPerPatch()
Returns the number of local boundary coefficients in each patch.
Array< OneD, int > m_bndCondCoeffsToLocalCoeffsMap
Integer map of bnd cond coeffs to local coefficients.
Definition: AssemblyMap.h:389
virtual const Array< OneD, const int > & v_GetExtraDirEdges()
int GetNumLocalCoeffs() const
Returns the total number of local coefficients.
std::string GetVariable()
Retrieves the variable string.
const Array< OneD, const unsigned int > & GetNumLocalIntCoeffsPerPatch()
Returns the number of local interior coefficients in each patch.
virtual void v_Assemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
bool m_signChange
Flag indicating if modes require sign reversal.
Definition: AssemblyMap.h:378
void Assemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
NekDouble GetLocalToGlobalBndSign(const int i) const
Retrieve the sign change of a given local boundary mode.
virtual std::shared_ptr< AssemblyMap > v_LinearSpaceMap(const ExpList &locexp, GlobalSysSolnType solnType)
Generate a linear space mapping from existing mapping.
void AssembleBnd(const NekVector< NekDouble > &loc, NekVector< NekDouble > &global, int offset) const
virtual void v_LocalToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, bool useComm) const
PatchMapSharedPtr m_patchMapFromPrevLevel
Mapping information for previous level in MultiLevel Solver.
Definition: AssemblyMap.h:514
Array< OneD, int > m_localToLocalIntMap
Integer map of local boundary coeffs to local interior system numbering.
Definition: AssemblyMap.h:387
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMapUnique()
NekDouble m_iterativeTolerance
Tolerance for iterative solver.
Definition: AssemblyMap.h:410
NekDouble GetIterativeTolerance() const
int m_numGlobalCoeffs
Total number of global coefficients.
Definition: AssemblyMap.h:375
std::string m_linSysIterSolver
Iterative solver: Conjugate Gradient, GMRES.
Definition: AssemblyMap.h:417
bool AtLastLevel() const
Returns true if this is the last level in the multi-level static condensation.
void UniversalAssembleBnd(Array< OneD, NekDouble > &pGlobal) const
const Array< OneD, const int > & GetLocalToGlobalMap()
const Array< OneD, const int > & GetBndCondCoeffsToLocalTraceMap()
Retrieves the local indices corresponding to the boundary expansion modes to global trace.
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMap()
Array< OneD, int > m_globalToUniversalBndMap
Integer map of process coeffs to universal space.
Definition: AssemblyMap.h:397
void PatchAssemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
int m_successiveRHS
sucessive RHS for iterative solver
Definition: AssemblyMap.h:414
std::shared_ptr< AssemblyMap > LinearSpaceMap(const ExpList &locexp, GlobalSysSolnType solnType)
virtual int v_GetNumNonDirVertexModes() const
void CalculateBndSystemBandWidth()
Calculates the bandwidth of the boundary system.
void LocalToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, bool useComm=true) const
Array< OneD, NekDouble > m_localToGlobalBndSign
Integer sign of local boundary coeffs to global space.
Definition: AssemblyMap.h:383
void UniversalAbsMaxBnd(Array< OneD, NekDouble > &bndvals)
Array< OneD, unsigned int > m_numLocalBndCoeffsPerPatch
The number of bnd dofs per patch.
Definition: AssemblyMap.h:430
virtual const Array< OneD, const int > & v_GetLocalToGlobalMap()
const Array< OneD, const int > & GetGlobalToUniversalMap()
std::string GetLinSysIterSolver() const
virtual int v_GetNumNonDirEdgeModes() const
virtual int v_GetNumDirFaces() const
LibUtilities::SessionReaderSharedPtr m_session
Session object.
Definition: AssemblyMap.h:333
virtual const Array< OneD, NekDouble > & v_GetLocalToGlobalSign() const
int GetNumGlobalDirBndCoeffs() const
Returns the number of global Dirichlet boundary coefficients.
const Array< OneD, NekDouble > & GetBndCondCoeffsToLocalCoeffsSign()
Returns the modal sign associated with a given boundary expansion mode.
int m_numLocalBndCoeffs
Number of local boundary coefficients.
Definition: AssemblyMap.h:345
void LocalToLocalBnd(const Array< OneD, const NekDouble > &local, Array< OneD, NekDouble > &locbnd) const
const Array< OneD, const int > & GetBndCondCoeffsToLocalCoeffsMap()
Retrieves the local indices corresponding to the boundary expansion modes.
void PrintStats(std::ostream &out, std::string variable, bool printHeader=true) const
AssemblyMapSharedPtr m_nextLevelLocalToGlobalMap
Map from the patches of the previous level to the patches of the current level.
Definition: AssemblyMap.h:437
const Array< OneD, const int > & GetLocalToGlobalBndMap()
Retrieve the global indices of the local boundary modes.
int m_staticCondLevel
The level of recursion in the case of multi-level static condensation.
Definition: AssemblyMap.h:426
int m_bndSystemBandWidth
The bandwith of the global bnd system.
Definition: AssemblyMap.h:404
const Array< OneD, const int > & GetGlobalToUniversalBndMapUnique()
void GlobalToLocal(const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const
GlobalSysSolnType GetGlobalSysSolnType() const
Returns the method of solving global systems.
int m_numLocalDirBndCoeffs
Number of Local Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:349
bool GetSingularSystem() const
Retrieves if the system is singular (true) or not (false)
std::string m_variable
Variable string identifier.
Definition: AssemblyMap.h:339
virtual int v_GetNumNonDirFaceModes() const
void PatchGlobalToLocal(const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const
virtual int v_GetFullSystemBandWidth() const
int m_numGlobalDirBndCoeffs
Number of Global Dirichlet Boundary Coefficients.
Definition: AssemblyMap.h:351
LibUtilities::CommSharedPtr GetComm()
Retrieves the communicator.
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
The number of int dofs per patch.
Definition: AssemblyMap.h:432
const AssemblyMapSharedPtr GetNextLevelLocalToGlobalMap() const
Returns the local to global mapping for the next level in the multi-level static condensation.
int GetBndSystemBandWidth() const
Returns the bandwidth of the boundary system.
virtual void v_UniversalAssemble(Array< OneD, NekDouble > &pGlobal) const
int GetStaticCondLevel() const
Returns the level of static condensation for this map.
void GlobalToLocalBndWithoutSign(const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc)
bool m_systemSingular
Flag indicating if the system is singular or not.
Definition: AssemblyMap.h:353
Array< OneD, int > m_localToGlobalBndMap
Integer map of local coeffs to global Boundary Dofs.
Definition: AssemblyMap.h:381
size_t GetHash() const
Retrieves the hash of this map.
int GetLocalToGlobalBndMap(const int i) const
Retrieve the global index of a given local boundary mode.
virtual ~AssemblyMap()
Destructor.
Gs::gs_data * m_dirBndGsh
gs gather communication to impose Dirhichlet BCs.
Definition: AssemblyMap.h:422
virtual void v_GlobalToLocal(const Array< OneD, const NekDouble > &global, Array< OneD, NekDouble > &loc) const
const Array< OneD, NekDouble > & GetLocalToGlobalSign() const
int GetNumLocalDirBndCoeffs() const
Returns the number of local Dirichlet boundary coefficients.
Array< OneD, const NekDouble > GetLocalToGlobalBndSign() const
Retrieve the sign change for all local boundary modes.
void LocalBndToLocal(const Array< OneD, const NekDouble > &locbnd, Array< OneD, NekDouble > &local) const
AssemblyMap()
Default constructor.
Definition: AssemblyMap.cpp:79
Array< OneD, int > m_globalToUniversalBndMapUnique
Integer map of unique process coeffs to universal space (signed)
Definition: AssemblyMap.h:399
virtual int v_GetNumNonDirEdges() const
virtual int v_GetNumDirEdges() const
const PatchMapSharedPtr & GetPatchMapFromPrevLevel(void) const
Returns the patch map from the previous level of the multi-level static condensation.
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: AssemblyMap.h:336
int GetNumGlobalBndCoeffs() const
Returns the total number of global boundary coefficients.
const Array< OneD, const int > & GetGlobalToUniversalBndMap()
Array< OneD, int > m_localToLocalBndMap
Integer map of local boundary coeffs to local boundary system numbering.
Definition: AssemblyMap.h:385
Array< OneD, NekDouble > m_bndCondCoeffsToLocalCoeffsSign
Integer map of sign of bnd cond coeffs to local coefficients.
Definition: AssemblyMap.h:391
const Array< OneD, const int > & GetBndCondIDToGlobalTraceID()
void LocalBndToGlobal(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global, int offset, bool UseComm=true) const
void UniversalAssemble(Array< OneD, NekDouble > &pGlobal) const
int m_numPatches
The number of patches (~elements) in the current level.
Definition: AssemblyMap.h:428
void LocalToLocalInt(const Array< OneD, const NekDouble > &local, Array< OneD, NekDouble > &locint) const
virtual int v_GetNumNonDirFaces() const
bool GetSignChange()
Returns true if using a modal expansion requiring a change of sign of some modes.
void LocalIntToLocal(const Array< OneD, const NekDouble > &locbnd, Array< OneD, NekDouble > &local) const
void GlobalToLocalBnd(const NekVector< NekDouble > &global, NekVector< NekDouble > &loc, int offset) const
std::string m_preconType
Type type of preconditioner to use in iterative solver.
Definition: AssemblyMap.h:407
int GetNumLocalBndCoeffs() const
Returns the total number of local boundary coefficients.
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:347
Array< OneD, int > m_bndCondIDToGlobalTraceID
Integer map of bnd cond trace number to global trace number.
Definition: AssemblyMap.h:395
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:98
Array< OneD, DataType > & GetPtr()
Definition: NekVector.cpp:217
static void Gather(Nektar::Array< OneD, NekDouble > pU, gs_op pOp, gs_data *pGsh, Nektar::Array< OneD, NekDouble > pBuffer=NullNekDouble1DArray)
Performs a gather-scatter operation of the provided values.
Definition: GsLib.hpp:265
@ gs_amax
Definition: GsLib.hpp:64
@ gs_add
Definition: GsLib.hpp:60
@ gs_max
Definition: GsLib.hpp:63
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
@ eNoSolnType
No Solution type specified.
std::shared_ptr< BottomUpSubStructuredGraph > BottomUpSubStructuredGraphSharedPtr
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:50
std::shared_ptr< PatchMap > PatchMapSharedPtr
int RoundNekDoubleToInt(NekDouble x)
Rounds a double precision number to an integer.
Definition: AssemblyMap.cpp:56
static const NekDouble kNekIterativeTol
double NekDouble
void Gathr(I n, const T *x, const I *y, T *z)
Gather vector z[i] = x[y[i]].
Definition: Vmath.hpp:507
void Scatr(int n, const T *x, const int *y, T *z)
Scatter vector z[y[i]] = x[i].
Definition: Vmath.hpp:539
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.hpp:577
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294