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, "AbsoluteTolerance"))
128 {
129 std::string abstol =
130 pSession->GetGlobalSysSolnInfo(variable, "AbsoluteTolerance");
132 boost::iequals(boost::to_upper_copy(abstol), "TRUE");
133 }
134 else
135 {
136 m_isAbsoluteTolerance = false;
137 }
138
139 if (pSession->DefinesGlobalSysSolnInfo(variable, "SuccessiveRHS"))
140 {
141 m_successiveRHS = boost::lexical_cast<int>(
142 pSession->GetGlobalSysSolnInfo(variable, "SuccessiveRHS").c_str());
143 }
144 else
145 {
146 pSession->LoadParameter("SuccessiveRHS", m_successiveRHS, 0);
147 }
148
149 if (pSession->DefinesGlobalSysSolnInfo(variable, "LinSysIterSolver"))
150 {
152 pSession->GetGlobalSysSolnInfo(variable, "LinSysIterSolver");
153 }
154 else if (pSession->DefinesSolverInfo("LinSysIterSolver"))
155 {
156 m_linSysIterSolver = pSession->GetSolverInfo("LinSysIterSolver");
157 }
158 else
159 {
160 m_linSysIterSolver = "ConjugateGradient";
161 }
162}
163
164/**
165 * Create a new level of mapping using the information in
166 * multiLevelGraph and performing the following steps:
167 */
169 AssemblyMap *oldLevelMap,
170 const BottomUpSubStructuredGraphSharedPtr &multiLevelGraph)
171 : m_session(oldLevelMap->m_session), m_comm(oldLevelMap->GetComm()),
172 m_hash(0), m_solnType(oldLevelMap->m_solnType),
173 m_preconType(oldLevelMap->m_preconType),
174 m_successiveRHS(oldLevelMap->m_successiveRHS),
175 m_linSysIterSolver(oldLevelMap->m_linSysIterSolver),
176 m_gsh(oldLevelMap->m_gsh), m_bndGsh(oldLevelMap->m_bndGsh),
177 m_lowestStaticCondLevel(oldLevelMap->m_lowestStaticCondLevel)
178{
179 int i;
180 int j;
181 int cnt;
182
183 //--------------------------------------------------------------
184 // -- Extract information from the input argument
185 int numGlobalBndCoeffsOld = oldLevelMap->GetNumGlobalBndCoeffs();
186 int numGlobalDirBndCoeffsOld = oldLevelMap->GetNumGlobalDirBndCoeffs();
187 int numLocalBndCoeffsOld = oldLevelMap->GetNumLocalBndCoeffs();
188 int numLocalDirBndCoeffsOld = oldLevelMap->GetNumLocalDirBndCoeffs();
189 bool signChangeOld = oldLevelMap->GetSignChange();
190
191 int staticCondLevelOld = oldLevelMap->GetStaticCondLevel();
192 int numPatchesOld = oldLevelMap->GetNumPatches();
193 GlobalSysSolnType solnTypeOld = oldLevelMap->GetGlobalSysSolnType();
194 const Array<OneD, const unsigned int> &numLocalBndCoeffsPerPatchOld =
195 oldLevelMap->GetNumLocalBndCoeffsPerPatch();
196 //--------------------------------------------------------------
197
198 //--------------------------------------------------------------
199 int newLevel = staticCondLevelOld + 1;
200 /** - STEP 1: setup a mask array to determine to which patch
201 * of the new level every patch of the current
202 * level belongs. To do so we make four arrays,
203 * #gloPatchMask, #globHomPatchMask,
204 * #locPatchMask_NekDouble and #locPatchMask.
205 * These arrays are then used to check which local
206 * dofs of the old level belong to which patch of
207 * the new level
208 */
209 Array<OneD, NekDouble> globPatchMask(numGlobalBndCoeffsOld, -1.0);
210 Array<OneD, NekDouble> globHomPatchMask(globPatchMask +
211 numGlobalDirBndCoeffsOld);
212 Array<OneD, NekDouble> locPatchMask_NekDouble(numLocalBndCoeffsOld, -3.0);
213 Array<OneD, int> locPatchMask(numLocalBndCoeffsOld);
214
215 // Fill the array globPatchMask as follows:
216 // - The first part (i.e. the glob bnd dofs) is filled with the
217 // value -1
218 // - The second part (i.e. the glob interior dofs) is numbered
219 // according to the patch it belongs to (i.e. dofs in first block
220 // all are numbered 0, the second block numbered are 1, etc...)
221 multiLevelGraph->MaskPatches(newLevel, globHomPatchMask);
222
223 // Map from Global Dofs to Local Dofs
224 // As a result, we know for each local dof whether
225 // it is mapped to the boundary of the next level, or to which
226 // patch. Based upon this, we can than later associate every patch
227 // of the current level with a patch in the next level.
228 oldLevelMap->GlobalToLocalBndWithoutSign(globPatchMask,
229 locPatchMask_NekDouble);
230
231 // Convert the result to an array of integers rather than doubles
232 RoundNekDoubleToInt(locPatchMask_NekDouble, locPatchMask);
233
234 /** - STEP 2: We calculate how many local bnd dofs of the
235 * old level belong to the boundaries of each patch at
236 * the new level. We need this information to set up the
237 * mapping between different levels.
238 */
239
240 // Retrieve the number of patches at the next level
241 int numPatchesWithIntNew =
242 multiLevelGraph->GetNpatchesWithInterior(newLevel);
243 int numPatchesNew = numPatchesWithIntNew;
244
245 // Allocate memory to store the number of local dofs associated to
246 // each of elemental boundaries of these patches
247 std::map<int, int> numLocalBndCoeffsPerPatchNew;
248 for (int i = 0; i < numPatchesNew; i++)
249 {
250 numLocalBndCoeffsPerPatchNew[i] = 0;
251 }
252
253 int minval;
254 int maxval;
255 int curPatch;
256 for (i = cnt = 0; i < numPatchesOld; i++)
257 {
258 // For every patch at the current level, the mask array
259 // locPatchMask should be filled with
260 // - the same (positive) number for each entry (which will
261 // correspond to the patch at the next level it belongs to)
262 // - the same (positive) number for each entry, except some
263 // entries that are -1 (the enties correspond to -1, will be
264 // mapped to the local boundary of the next level patch given
265 // by the positive number)
266 // - -1 for all entries. In this case, we will make an
267 // additional patch only consisting of boundaries at the next
268 // level
269 minval =
270 *min_element(&locPatchMask[cnt],
271 &locPatchMask[cnt] + numLocalBndCoeffsPerPatchOld[i]);
272 maxval =
273 *max_element(&locPatchMask[cnt],
274 &locPatchMask[cnt] + numLocalBndCoeffsPerPatchOld[i]);
275 ASSERTL0((minval == maxval) || (minval == -1),
276 "These values should never be the same");
277
278 if (maxval == -1)
279 {
280 curPatch = numPatchesNew;
281 numLocalBndCoeffsPerPatchNew[curPatch] = 0;
282 numPatchesNew++;
283 }
284 else
285 {
286 curPatch = maxval;
287 }
288
289 for (j = 0; j < numLocalBndCoeffsPerPatchOld[i]; j++)
290 {
291 ASSERTL0((locPatchMask[cnt] == maxval) ||
292 (locPatchMask[cnt] == minval),
293 "These values should never be the same");
294 if (locPatchMask[cnt] == -1)
295 {
296 ++numLocalBndCoeffsPerPatchNew[curPatch];
297 }
298 cnt++;
299 }
300 }
301
302 // Count how many local dofs of the old level are mapped
303 // to the local boundary dofs of the new level
306 m_numPatches = numLocalBndCoeffsPerPatchNew.size();
309 multiLevelGraph->GetNintDofsPerPatch(newLevel, m_numLocalIntCoeffsPerPatch);
310
311 for (int i = 0; i < m_numPatches; i++)
312 {
314 (unsigned int)numLocalBndCoeffsPerPatchNew[i];
318 }
319
320 // Also initialise some more data members
321 m_solnType = solnTypeOld;
326 "This method should only be called for in "
327 "case of multi-level static condensation.");
328 m_staticCondLevel = newLevel;
329 m_signChange = signChangeOld;
330 m_numLocalDirBndCoeffs = numLocalDirBndCoeffsOld;
331 m_numGlobalDirBndCoeffs = numGlobalDirBndCoeffsOld;
333 multiLevelGraph->GetInteriorOffset(newLevel) + m_numGlobalDirBndCoeffs;
335 multiLevelGraph->GetNumGlobalDofs(newLevel) + m_numGlobalDirBndCoeffs;
337
341
342 // local to bnd map is just a copy
343 for (int i = 0; i < m_numLocalBndCoeffs; ++i)
344 {
346 }
347
348 // local to int map is just a copy plus offset
349 for (int i = m_numLocalBndCoeffs; i < m_numLocalCoeffs; ++i)
350 {
352 }
353
354 if (m_signChange)
355 {
357 }
358
360 MemoryManager<PatchMap>::AllocateSharedPtr(numLocalBndCoeffsOld);
361
366
367 // Set up an offset array that denotes the offset of the local
368 // boundary degrees of freedom of the next level
369 Array<OneD, int> numLocalBndCoeffsPerPatchOffset(m_numPatches + 1, 0);
370 for (int i = 1; i < m_numPatches + 1; i++)
371 {
372 numLocalBndCoeffsPerPatchOffset[i] +=
373 numLocalBndCoeffsPerPatchOffset[i - 1] +
374 numLocalBndCoeffsPerPatchNew[i - 1];
375 }
376
377 int additionalPatchCnt = numPatchesWithIntNew;
378 int newid;
379 int blockid;
380 bool isBndDof;
382 Array<OneD, int> bndDofPerPatchCnt(m_numPatches, 0);
383
384 for (i = cnt = 0; i < numPatchesOld; i++)
385 {
386 minval =
387 *min_element(&locPatchMask[cnt],
388 &locPatchMask[cnt] + numLocalBndCoeffsPerPatchOld[i]);
389 maxval =
390 *max_element(&locPatchMask[cnt],
391 &locPatchMask[cnt] + numLocalBndCoeffsPerPatchOld[i]);
392 ASSERTL0((minval == maxval) || (minval == -1),
393 "These values should never be the same");
394
395 if (maxval == -1)
396 {
397 curPatch = additionalPatchCnt;
398 additionalPatchCnt++;
399 }
400 else
401 {
402 curPatch = maxval;
403 }
404
405 for (j = 0; j < numLocalBndCoeffsPerPatchOld[i]; j++)
406 {
407 ASSERTL0((locPatchMask[cnt] == maxval) ||
408 (locPatchMask[cnt] == minval),
409 "These values should never be the same");
410
411 sign = oldLevelMap->GetLocalToGlobalBndSign(cnt);
412
413 if (locPatchMask[cnt] == -1)
414 {
415 newid = numLocalBndCoeffsPerPatchOffset[curPatch];
416
417 m_localToGlobalBndMap[newid] =
418 oldLevelMap->GetLocalToGlobalBndMap(cnt);
419
420 if (m_signChange)
421 {
423 }
424
425 blockid = bndDofPerPatchCnt[curPatch];
426 isBndDof = true;
427
428 numLocalBndCoeffsPerPatchOffset[curPatch]++;
429 bndDofPerPatchCnt[curPatch]++;
430 }
431 else
432 {
433 newid = oldLevelMap->GetLocalToGlobalBndMap(cnt) -
435
436 blockid =
437 oldLevelMap->GetLocalToGlobalBndMap(cnt) -
439 multiLevelGraph->GetInteriorOffset(newLevel, curPatch);
440 isBndDof = false;
441 }
442
443 sign = isBndDof ? 1.0 : sign;
444
445 m_patchMapFromPrevLevel->SetPatchMap(cnt, curPatch, blockid,
446 isBndDof, sign);
447 cnt++;
448 }
449 }
450
451 // set up local to local mapping from previous to new level
454
456
457 // Postprocess the computed information - Update the old
458 // level with the mapping to new level
459 // oldLevelMap->SetLocalBndToNextLevelMap(oldLocalBndToNewLevelMap,oldLocalBndToNewLevelSign);
460
461 // - Construct the next level mapping object
462 if (m_staticCondLevel < (multiLevelGraph->GetNlevels() - 1))
463 {
466 multiLevelGraph);
467 }
468}
469
471{
472}
473
474/**
475 * The bandwidth calculated corresponds to what is referred to as
476 * half-bandwidth. If the elements of the matrix are designated as
477 * a_ij, it corresponds to the maximum value of |i-j| for non-zero
478 * a_ij. As a result, the value also corresponds to the number of
479 * sub- or super-diagonals.
480 *
481 * The bandwith can be calculated elementally as it corresponds to the
482 * maximal elemental bandwith (i.e. the maximal difference in global
483 * DOF index for every element).
484 *
485 * We here calculate the bandwith of the global boundary system (as
486 * used for static condensation).
487 */
489{
490 int i, j;
491 int cnt = 0;
492 int locSize;
493 int maxId;
494 int minId;
495 int bwidth = -1;
496 for (i = 0; i < m_numPatches; ++i)
497 {
498 locSize = m_numLocalBndCoeffsPerPatch[i];
499 maxId = -1;
500 minId = m_numLocalBndCoeffs + 1;
501 for (j = 0; j < locSize; j++)
502 {
504 {
505 if (m_localToGlobalBndMap[cnt + j] > maxId)
506 {
507 maxId = m_localToGlobalBndMap[cnt + j];
508 }
509
510 if (m_localToGlobalBndMap[cnt + j] < minId)
511 {
512 minId = m_localToGlobalBndMap[cnt + j];
513 }
514 }
515 }
516 bwidth = (bwidth > (maxId - minId)) ? bwidth : (maxId - minId);
517
518 cnt += locSize;
519 }
520
521 m_bndSystemBandWidth = bwidth;
522}
523
524int AssemblyMap::v_GetLocalToGlobalMap([[maybe_unused]] const int i) const
525{
526 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
527 return 0;
528}
529
530int AssemblyMap::v_GetGlobalToUniversalMap([[maybe_unused]] const int i) const
531{
532 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
533 return 0;
534}
535
537 [[maybe_unused]] const int i) const
538{
539 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
540 return 0;
541}
542
544{
545 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
546 static Array<OneD, const int> result;
547 return result;
548}
549
551{
552 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
553 static Array<OneD, const int> result;
554 return result;
555}
556
558{
559 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
560 static Array<OneD, const int> result;
561 return result;
562}
563
565 [[maybe_unused]] const int i) const
566{
567 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
568 return 0.0;
569}
570
572{
573 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
574 static Array<OneD, NekDouble> result;
575 return result;
576}
577
579 [[maybe_unused]] const Array<OneD, const NekDouble> &loc,
580 [[maybe_unused]] Array<OneD, NekDouble> &global,
581 [[maybe_unused]] bool useComm) const
582{
583 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
584}
585
587 [[maybe_unused]] const NekVector<NekDouble> &loc,
588 [[maybe_unused]] NekVector<NekDouble> &global,
589 [[maybe_unused]] bool useComm) const
590{
591 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
592}
593
595 [[maybe_unused]] const Array<OneD, const NekDouble> &global,
596 [[maybe_unused]] Array<OneD, NekDouble> &loc) const
597{
598 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
599}
600
602 [[maybe_unused]] const NekVector<NekDouble> &global,
603 [[maybe_unused]] NekVector<NekDouble> &loc) const
604{
605 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
606}
607
609 [[maybe_unused]] const Array<OneD, const NekDouble> &loc,
610 [[maybe_unused]] Array<OneD, NekDouble> &global) const
611{
612 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
613}
614
616 [[maybe_unused]] const NekVector<NekDouble> &loc,
617 [[maybe_unused]] NekVector<NekDouble> &global) const
618{
619 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
620}
621
623 [[maybe_unused]] Array<OneD, NekDouble> &pGlobal) const
624{
625 // Do nothing here since multi-level static condensation uses a
626 // AssemblyMap and thus will call this routine in serial.
627}
628
630 [[maybe_unused]] NekVector<NekDouble> &pGlobal) const
631{
632 // Do nothing here since multi-level static condensation uses a
633 // AssemblyMap and thus will call this routine in serial.
634}
635
637 [[maybe_unused]] Array<OneD, NekDouble> &pGlobal,
638 [[maybe_unused]] int offset) 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{
646 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
647 return 0;
648}
649
651{
652 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
653 return 0;
654}
655
657{
658 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
659 return 0;
660}
661
663{
664 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
665 return 0;
666}
667
669{
670 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
671 return 0;
672}
673
675{
676 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
677 return 0;
678}
679
681{
682 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
683 return 0;
684}
685
687{
688 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
689 return 0;
690}
691
693{
694 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
695 static Array<OneD, const int> result;
696 return result;
697}
698
699std::shared_ptr<AssemblyMap> AssemblyMap::v_LinearSpaceMap(
700 [[maybe_unused]] const ExpList &locexp,
701 [[maybe_unused]] GlobalSysSolnType solnType)
702{
703 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
704 static std::shared_ptr<AssemblyMap> result;
705 return result;
706}
707
709{
710 return m_comm;
711}
712
714{
715 return m_variable;
716}
717
719{
720 return m_hash;
721}
722
724{
725 return v_GetLocalToGlobalMap(i);
726}
727
729{
731}
732
734{
736}
737
739{
740 return v_GetLocalToGlobalMap();
741}
742
744{
746}
747
749{
751}
752
754{
755 return v_GetLocalToGlobalSign(i);
756}
757
759{
760 return v_GetLocalToGlobalSign();
761}
762
765 bool useComm) const
766{
767 v_LocalToGlobal(loc, global, useComm);
768}
769
771 NekVector<NekDouble> &global,
772 bool useComm) const
773{
774 v_LocalToGlobal(loc, global, useComm);
775}
776
779{
780 v_GlobalToLocal(global, loc);
781}
782
785{
786 v_GlobalToLocal(global, loc);
787}
788
790 Array<OneD, NekDouble> &global) const
791{
792 v_Assemble(loc, global);
793}
794
796 NekVector<NekDouble> &global) const
797{
798 v_Assemble(loc, global);
799}
800
802{
803 v_UniversalAssemble(pGlobal);
804}
805
807{
808 v_UniversalAssemble(pGlobal);
809}
810
812 int offset) const
813{
814 v_UniversalAssemble(pGlobal, offset);
815}
816
818 Array<OneD, NekDouble> &global) const
819{
821
822 Array<OneD, const int> map = m_patchMapFromPrevLevel->GetNewLevelMap();
824
825 if (global.data() == loc.data())
826 {
827 local = Array<OneD, NekDouble>(map.size(), loc.data());
828 }
829 else
830 {
831 local = loc; // create reference
832 }
833
834 Vmath::Scatr(map.size(), sign.get(), local.get(), map.get(), global.get());
835}
836
839{
841
842 Array<OneD, const int> map = m_patchMapFromPrevLevel->GetNewLevelMap();
844
845 if (global.data() == loc.data())
846 {
847 glo = Array<OneD, NekDouble>(global.size(), global.data());
848 }
849 else
850 {
851 glo = global; // create reference
852 }
853
854 Vmath::Gathr(map.size(), sign.get(), glo.get(), map.get(), loc.get());
855}
856
858 Array<OneD, NekDouble> &global) const
859{
861 Array<OneD, const int> map = m_patchMapFromPrevLevel->GetNewLevelMap();
863
864 if (global.data() == loc.data())
865 {
866 local = Array<OneD, NekDouble>(map.size(), loc.data());
867 }
868 else
869 {
870 local = loc; // create reference
871 }
872
873 // since we are calling mapping from level down from array
874 // the m_numLocaBndCoeffs represents the size of the
875 // boundary elements we need to assemble into
876 Vmath::Zero(m_numLocalCoeffs, global.get(), 1);
877
878 Vmath::Assmb(map.size(), sign.get(), local.get(), map.get(), global.get());
879}
880
882{
884}
885
887{
889}
890
892{
894}
895
897{
899}
900
902{
903 return v_GetNumDirEdges();
904}
905
907{
908 return v_GetNumDirFaces();
909}
910
912{
913 return v_GetNumNonDirEdges();
914}
915
917{
918 return v_GetNumNonDirFaces();
919}
920
922{
923 return v_GetExtraDirEdges();
924}
925
926std::shared_ptr<AssemblyMap> AssemblyMap::LinearSpaceMap(
927 const ExpList &locexp, GlobalSysSolnType solnType)
928{
929 return v_LinearSpaceMap(locexp, solnType);
930}
931
933{
934 return m_localToGlobalBndMap[i];
935}
936
938{
940}
941
943{
944 return m_signChange;
945}
946
948{
950}
951
953{
955}
956
958{
960}
961
963{
964 if (m_signChange)
965 {
966 return m_localToGlobalBndSign[i];
967 }
968 else
969 {
970 return 1.0;
971 }
972}
973
975{
977}
978
980{
982}
983
985{
987}
988
990{
991 ASSERTL1(i < m_bndCondIDToGlobalTraceID.size(), "Index out of range.");
993}
994
996{
998}
999
1001{
1003}
1004
1006{
1008}
1009
1011{
1012 return m_numLocalBndCoeffs;
1013}
1014
1016{
1017 return m_numGlobalBndCoeffs;
1018}
1019
1021{
1022 return m_numLocalCoeffs;
1023}
1024
1026{
1027 return m_numGlobalCoeffs;
1028}
1029
1031{
1032 return m_systemSingular;
1033}
1034
1036 NekVector<NekDouble> &loc, int offset) const
1037{
1038 GlobalToLocalBnd(global.GetPtr(), loc.GetPtr(), offset);
1039}
1040
1043{
1044 GlobalToLocalBnd(global.GetPtr(), loc.GetPtr());
1045}
1046
1049 int offset) const
1050{
1052 "Local vector is not of correct dimension");
1053 ASSERTL1(global.size() >= m_numGlobalBndCoeffs - offset,
1054 "Global vector is not of correct dimension");
1055
1056 // offset input data by length "offset" for Dirichlet boundary conditions.
1058 Vmath::Vcopy(m_numGlobalBndCoeffs - offset, global.get(), 1,
1059 tmp.get() + offset, 1);
1060
1061 if (m_signChange)
1062 {
1064 tmp.get(), m_localToGlobalBndMap.get(), loc.get());
1065 }
1066 else
1067 {
1069 m_localToGlobalBndMap.get(), loc.get());
1070 }
1071}
1072
1075{
1077 "Local vector is not of correct dimension");
1078 ASSERTL1(global.size() >= m_numGlobalBndCoeffs,
1079 "Global vector is not of correct dimension");
1080
1082 if (global.data() == loc.data())
1083 {
1084 glo = Array<OneD, NekDouble>(m_numLocalBndCoeffs, global.data());
1085 }
1086 else
1087 {
1088 glo = global; // create reference
1089 }
1090
1091 if (m_signChange)
1092 {
1094 glo.get(), m_localToGlobalBndMap.get(), loc.get());
1095 }
1096 else
1097 {
1099 m_localToGlobalBndMap.get(), loc.get());
1100 }
1101}
1102
1104 Array<OneD, NekDouble> &global, int offset,
1105 bool UseComm) const
1106{
1108 "Local vector is not of correct dimension");
1109 ASSERTL1(global.size() >= m_numGlobalBndCoeffs - offset,
1110 "Global vector is not of correct dimension");
1111
1112 // offset input data by length "offset" for Dirichlet boundary conditions.
1114
1115 if (m_signChange)
1116 {
1118 loc.get(), m_localToGlobalBndMap.get(), tmp.get());
1119 }
1120 else
1121 {
1123 m_localToGlobalBndMap.get(), tmp.get());
1124 }
1125
1126 // Ensure each processor has unique value with a max gather.
1127 if (UseComm)
1128 {
1130 }
1131 Vmath::Vcopy(m_numGlobalBndCoeffs - offset, tmp.get() + offset, 1,
1132 global.get(), 1);
1133}
1134
1136 Array<OneD, NekDouble> &global,
1137 bool UseComm) const
1138{
1140 "Local vector is not of correct dimension");
1141 ASSERTL1(global.size() >= m_numGlobalBndCoeffs,
1142 "Global vector is not of correct dimension");
1143
1144 if (m_signChange)
1145 {
1147 loc.get(), m_localToGlobalBndMap.get(), global.get());
1148 }
1149 else
1150 {
1152 m_localToGlobalBndMap.get(), global.get());
1153 }
1154 if (UseComm)
1155 {
1156 Gs::Gather(global, Gs::gs_max, m_bndGsh);
1157 }
1158}
1159
1161 Array<OneD, NekDouble> &locbnd) const
1162{
1163 ASSERTL1(locbnd.size() >= m_numLocalBndCoeffs,
1164 "LocBnd vector is not of correct dimension");
1165 ASSERTL1(local.size() >= m_numLocalCoeffs,
1166 "Local vector is not of correct dimension");
1167
1169 locbnd.get());
1170}
1171
1173 Array<OneD, NekDouble> &locint) const
1174{
1176 "Locint vector is not of correct dimension");
1177 ASSERTL1(local.size() >= m_numLocalCoeffs,
1178 "Local vector is not of correct dimension");
1179
1181 m_localToLocalIntMap.get(), locint.get());
1182}
1183
1185 Array<OneD, NekDouble> &local) const
1186{
1187 ASSERTL1(locbnd.size() >= m_numLocalBndCoeffs,
1188 "LocBnd vector is not of correct dimension");
1189 ASSERTL1(local.size() >= m_numLocalCoeffs,
1190 "Local vector is not of correct dimension");
1191
1193 local.get());
1194}
1195
1197 Array<OneD, NekDouble> &local) const
1198{
1200 "LocBnd vector is not of correct dimension");
1201 ASSERTL1(local.size() >= m_numLocalCoeffs,
1202 "Local vector is not of correct dimension");
1203
1205 m_localToLocalIntMap.get(), local.get());
1206}
1207
1209 NekVector<NekDouble> &global, int offset) const
1210{
1211 AssembleBnd(loc.GetPtr(), global.GetPtr(), offset);
1212}
1213
1215 NekVector<NekDouble> &global) const
1216{
1217 AssembleBnd(loc.GetPtr(), global.GetPtr());
1218}
1219
1221 Array<OneD, NekDouble> &global, int offset) const
1222{
1224 "Local array is not of correct dimension");
1225 ASSERTL1(global.size() >= m_numGlobalBndCoeffs - offset,
1226 "Global array is not of correct dimension");
1228
1229 if (m_signChange)
1230 {
1232 loc.get(), m_localToGlobalBndMap.get(), tmp.get());
1233 }
1234 else
1235 {
1237 m_localToGlobalBndMap.get(), tmp.get());
1238 }
1240 Vmath::Vcopy(m_numGlobalBndCoeffs - offset, tmp.get() + offset, 1,
1241 global.get(), 1);
1242}
1243
1245 Array<OneD, NekDouble> &global) const
1246{
1248 "Local vector is not of correct dimension");
1249 ASSERTL1(global.size() >= m_numGlobalBndCoeffs,
1250 "Global vector is not of correct dimension");
1251
1252 Vmath::Zero(m_numGlobalBndCoeffs, global.get(), 1);
1253
1254 if (m_signChange)
1255 {
1257 loc.get(), m_localToGlobalBndMap.get(), global.get());
1258 }
1259 else
1260 {
1262 m_localToGlobalBndMap.get(), global.get());
1263 }
1264 UniversalAssembleBnd(global);
1265}
1266
1268{
1269 ASSERTL1(pGlobal.size() >= m_numGlobalBndCoeffs, "Wrong size.");
1270 Gs::Gather(pGlobal, Gs::gs_add, m_bndGsh);
1271}
1272
1274{
1275 UniversalAssembleBnd(pGlobal.GetPtr());
1276}
1277
1279 int offset) const
1280{
1281 Array<OneD, NekDouble> tmp(offset);
1282 if (offset > 0)
1283 {
1284 Vmath::Vcopy(offset, pGlobal, 1, tmp, 1);
1285 }
1286 UniversalAssembleBnd(pGlobal);
1287 if (offset > 0)
1288 {
1289 Vmath::Vcopy(offset, tmp, 1, pGlobal, 1);
1290 }
1291}
1292
1294{
1296}
1297
1299{
1300 return m_bndSystemBandWidth;
1301}
1302
1304{
1305 return m_staticCondLevel;
1306}
1307
1309{
1310 return m_numPatches;
1311}
1312
1315{
1317}
1318
1321{
1323}
1324
1326{
1328}
1329
1331{
1333}
1334
1336{
1337 return !(m_nextLevelLocalToGlobalMap.get());
1338}
1339
1341{
1342 return m_solnType;
1343}
1344
1346{
1347 return m_preconType;
1348}
1349
1351{
1352 return m_isAbsoluteTolerance;
1353}
1354
1356{
1357 return m_successiveRHS;
1358}
1359
1361{
1362 return m_linSysIterSolver;
1363}
1364
1367{
1369 "Local vector is not of correct dimension");
1370 ASSERTL1(global.size() >= m_numGlobalBndCoeffs,
1371 "Global vector is not of correct dimension");
1372
1374 loc.get());
1375}
1376
1377void AssemblyMap::PrintStats(std::ostream &out, std::string variable,
1378 bool printHeader) const
1379{
1380 LibUtilities::CommSharedPtr vRowComm = m_session->GetComm()->GetRowComm();
1381 bool isRoot = m_session->GetComm()->IsParallelInTime()
1382 ? m_session->GetComm()->GetRank() == 0
1383 : vRowComm->GetRank() == 0;
1384 int n = vRowComm->GetSize();
1385 int i;
1386
1387 // Determine number of global degrees of freedom.
1388 int globBndCnt = 0, globDirCnt = 0;
1389
1390 for (i = 0; i < m_numGlobalBndCoeffs; ++i)
1391 {
1393 {
1394 globBndCnt++;
1395
1397 {
1398 globDirCnt++;
1399 }
1400 }
1401 }
1402
1403 int globCnt = m_numGlobalCoeffs - m_numGlobalBndCoeffs + globBndCnt;
1404
1405 // Calculate maximum valency
1408
1410 tmpGlob.get());
1411 UniversalAssembleBnd(tmpGlob);
1412
1413 int totGlobDof = globCnt;
1414 int totGlobBndDof = globBndCnt;
1415 int totGlobDirDof = globDirCnt;
1416 int totLocalDof = m_numLocalCoeffs;
1417 int totLocalBndDof = m_numLocalBndCoeffs;
1418 int totLocalDirDof = m_numLocalDirBndCoeffs;
1419
1420 int meanValence = 0;
1421 int maxValence = 0;
1422 int minValence = 10000000;
1423 for (int i = 0; i < m_numGlobalBndCoeffs; ++i)
1424 {
1426 {
1427 continue;
1428 }
1429
1430 if (tmpGlob[i] > maxValence)
1431 {
1432 maxValence = tmpGlob[i];
1433 }
1434 if (tmpGlob[i] < minValence)
1435 {
1436 minValence = tmpGlob[i];
1437 }
1438 meanValence += tmpGlob[i];
1439 }
1440
1441 vRowComm->AllReduce(maxValence, LibUtilities::ReduceMax);
1442 vRowComm->AllReduce(minValence, LibUtilities::ReduceMin);
1443 vRowComm->AllReduce(meanValence, LibUtilities::ReduceSum);
1444 vRowComm->AllReduce(totGlobDof, LibUtilities::ReduceSum);
1445 vRowComm->AllReduce(totGlobBndDof, LibUtilities::ReduceSum);
1446 vRowComm->AllReduce(totGlobDirDof, LibUtilities::ReduceSum);
1447 vRowComm->AllReduce(totLocalDof, LibUtilities::ReduceSum);
1448 vRowComm->AllReduce(totLocalBndDof, LibUtilities::ReduceSum);
1449 vRowComm->AllReduce(totLocalDirDof, LibUtilities::ReduceSum);
1450
1451 meanValence /= totGlobBndDof;
1452
1453 if (isRoot)
1454 {
1455 if (printHeader)
1456 {
1457 out << "Assembly map statistics for field " << variable << ":"
1458 << endl;
1459 }
1460
1461 out << " - Number of local/global dof : " << totLocalDof
1462 << " " << totGlobDof << endl;
1463 out << " - Number of local/global boundary dof : " << totLocalBndDof
1464 << " " << totGlobBndDof << endl;
1465 out << " - Number of local/global Dirichlet dof : " << totLocalDirDof
1466 << " " << totGlobDirDof << endl;
1467 out << " - dof valency (min/max/mean) : " << minValence
1468 << " " << maxValence << " " << meanValence << endl;
1469
1470 if (n > 1)
1471 {
1472 NekDouble mean = m_numLocalCoeffs, mean2 = mean * mean;
1473 NekDouble minval = mean, maxval = mean;
1475
1476 for (i = 1; i < n; ++i)
1477 {
1478 vRowComm->Recv(i, tmp);
1479 mean += tmp[0];
1480 mean2 += tmp[0] * tmp[0];
1481
1482 if (tmp[0] > maxval)
1483 {
1484 maxval = tmp[0];
1485 }
1486 if (tmp[0] < minval)
1487 {
1488 minval = tmp[0];
1489 }
1490 }
1491
1492 if (maxval > 0.1)
1493 {
1494 out << " - Local dof dist. (min/max/mean/dev) : " << minval
1495 << " " << maxval << " " << (mean / n) << " "
1496 << sqrt(mean2 / n - mean * mean / n / n) << endl;
1497 }
1498
1499 vRowComm->Block();
1500
1501 mean = minval = maxval = m_numLocalBndCoeffs;
1502 mean2 = mean * mean;
1503
1504 for (i = 1; i < n; ++i)
1505 {
1506 vRowComm->Recv(i, tmp);
1507 mean += tmp[0];
1508 mean2 += tmp[0] * tmp[0];
1509
1510 if (tmp[0] > maxval)
1511 {
1512 maxval = tmp[0];
1513 }
1514 if (tmp[0] < minval)
1515 {
1516 minval = tmp[0];
1517 }
1518 }
1519
1520 out << " - Local bnd dof dist. (min/max/mean/dev) : " << minval
1521 << " " << maxval << " " << (mean / n) << " "
1522 << sqrt(mean2 / n - mean * mean / n / n) << endl;
1523 }
1524 }
1525 else
1526 {
1528 tmp[0] = m_numLocalCoeffs;
1529 vRowComm->Send(0, tmp);
1530 vRowComm->Block();
1531 tmp[0] = m_numLocalBndCoeffs;
1532 vRowComm->Send(0, tmp);
1533 }
1534
1535 // Either we have no more levels in the static condensation, or we
1536 // are not multi-level.
1538 {
1539 return;
1540 }
1541
1542 int level = 2;
1544 while (tmp->m_nextLevelLocalToGlobalMap)
1545 {
1546 tmp = tmp->m_nextLevelLocalToGlobalMap;
1547 ++level;
1548 }
1549
1550 // Print out multi-level static condensation information.
1551 if (n > 1)
1552 {
1553 if (isRoot)
1554 {
1555 NekDouble mean = level, mean2 = mean * mean;
1556 int minval = level, maxval = level;
1557
1558 Array<OneD, NekDouble> tmpRecv(1);
1559 for (i = 1; i < n; ++i)
1560 {
1561 vRowComm->Recv(i, tmpRecv);
1562 mean += tmpRecv[0];
1563 mean2 += tmpRecv[0] * tmpRecv[0];
1564
1565 if (tmpRecv[0] > maxval)
1566 {
1567 maxval = (int)(tmpRecv[0] + 0.5);
1568 }
1569 if (tmpRecv[0] < minval)
1570 {
1571 minval = (int)(tmpRecv[0] + 0.5);
1572 }
1573 }
1574
1575 out << " - M-level sc. dist. (min/max/mean/dev) : " << minval
1576 << " " << maxval << " " << (mean / n) << " "
1577 << sqrt(mean2 / n - mean * mean / n / n) << endl;
1578 }
1579 else
1580 {
1581 Array<OneD, NekDouble> tmpSend(1);
1582 tmpSend[0] = level;
1583 vRowComm->Send(0, tmpSend);
1584 }
1585 }
1586 else
1587 {
1588 if (isRoot)
1589 {
1590 out << " - Number of static cond. levels : " << level
1591 << endl;
1592 }
1593 }
1594
1595 if (isRoot)
1596 {
1597 out << "Stats at lowest static cond. level:" << endl;
1598 }
1599 tmp->PrintStats(out, variable, false);
1600}
1601} // 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:399
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:361
Array< OneD, int > m_bndCondCoeffsToLocalTraceMap
Integer map of bnd cond coeff to local trace coeff.
Definition: AssemblyMap.h:390
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:386
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:375
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:511
Array< OneD, int > m_localToLocalIntMap
Integer map of local boundary coeffs to local interior system numbering.
Definition: AssemblyMap.h:384
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMapUnique()
int m_numGlobalCoeffs
Total number of global coefficients.
Definition: AssemblyMap.h:372
std::string m_linSysIterSolver
Iterative solver: Conjugate Gradient, GMRES.
Definition: AssemblyMap.h:414
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:394
void PatchAssemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
int m_successiveRHS
sucessive RHS for iterative solver
Definition: AssemblyMap.h:411
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:380
void UniversalAbsMaxBnd(Array< OneD, NekDouble > &bndvals)
Array< OneD, unsigned int > m_numLocalBndCoeffsPerPatch
The number of bnd dofs per patch.
Definition: AssemblyMap.h:427
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:330
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:342
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:434
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:423
int m_bndSystemBandWidth
The bandwith of the global bnd system.
Definition: AssemblyMap.h:401
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:346
bool GetSingularSystem() const
Retrieves if the system is singular (true) or not (false)
std::string m_variable
Variable string identifier.
Definition: AssemblyMap.h:336
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:348
LibUtilities::CommSharedPtr GetComm()
Retrieves the communicator.
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
The number of int dofs per patch.
Definition: AssemblyMap.h:429
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:350
Array< OneD, int > m_localToGlobalBndMap
Integer map of local coeffs to global Boundary Dofs.
Definition: AssemblyMap.h:378
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:419
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:396
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:333
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:382
Array< OneD, NekDouble > m_bndCondCoeffsToLocalCoeffsSign
Integer map of sign of bnd cond coeffs to local coefficients.
Definition: AssemblyMap.h:388
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:425
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:404
int GetNumLocalBndCoeffs() const
Returns the total number of local boundary coefficients.
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:344
Array< OneD, int > m_bndCondIDToGlobalTraceID
Integer map of bnd cond trace number to global trace number.
Definition: AssemblyMap.h:392
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:99
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:278
@ 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
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
STL namespace.
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294