Nektar++
Loading...
Searching...
No Matches
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 = std::stoi(
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
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 Array<OneD, const NekDouble> &global,
588 [[maybe_unused]] Array<OneD, NekDouble> &loc) const
589{
590 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
591}
592
594 [[maybe_unused]] const NekVector<NekDouble> &global,
595 [[maybe_unused]] NekVector<NekDouble> &loc) const
596{
597 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
598}
599
601 [[maybe_unused]] const Array<OneD, const NekDouble> &loc,
602 [[maybe_unused]] Array<OneD, NekDouble> &global) const
603{
604 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
605}
606
608 [[maybe_unused]] const NekVector<NekDouble> &loc,
609 [[maybe_unused]] NekVector<NekDouble> &global) const
610{
611 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
612}
613
615 [[maybe_unused]] Array<OneD, NekDouble> &pGlobal) const
616{
617 // Do nothing here since multi-level static condensation uses a
618 // AssemblyMap and thus will call this routine in serial.
619}
620
622 [[maybe_unused]] Array<OneD, NekDouble> &pGlobal,
623 [[maybe_unused]] int offset) 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{
631 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
632 return 0;
633}
634
636{
637 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
638 return 0;
639}
640
642{
643 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
644 return 0;
645}
646
648{
649 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
650 return 0;
651}
652
654{
655 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
656 return 0;
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 static Array<OneD, const int> result;
681 return result;
682}
683
684std::shared_ptr<AssemblyMap> AssemblyMap::v_LinearSpaceMap(
685 [[maybe_unused]] const ExpList &locexp,
686 [[maybe_unused]] GlobalSysSolnType solnType)
687{
688 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
689 static std::shared_ptr<AssemblyMap> result;
690 return result;
691}
692
697
699{
700 return m_variable;
701}
702
704{
705 return m_hash;
706}
707
709{
710 return v_GetLocalToGlobalMap(i);
711}
712
714{
716}
717
722
727
732
737
739{
740 return v_GetLocalToGlobalSign(i);
741}
742
747
750 bool useComm) const
751{
752 v_LocalToGlobal(loc, global, useComm);
753}
754
756 NekVector<NekDouble> &global,
757 bool useComm) const
758{
759 v_LocalToGlobal(loc.GetPtr(), global.GetPtr(), useComm);
760}
761
763 Array<OneD, NekDouble> &loc) const
764{
765 v_GlobalToLocal(global, loc);
766}
767
769 NekVector<NekDouble> &loc) const
770{
771 v_GlobalToLocal(global, loc);
772}
773
775 Array<OneD, NekDouble> &global) const
776{
777 v_Assemble(loc, global);
778}
779
781 NekVector<NekDouble> &global) const
782{
783 v_Assemble(loc, global);
784}
785
790
795
797 int offset) const
798{
799 v_UniversalAssemble(pGlobal, offset);
800}
801
803 Array<OneD, NekDouble> &global) const
804{
806
807 Array<OneD, const int> map = m_patchMapFromPrevLevel->GetNewLevelMap();
809
810 if (global.data() == loc.data())
811 {
812 local = Array<OneD, NekDouble>(map.size(), loc.data());
813 }
814 else
815 {
816 local = loc; // create reference
817 }
818
819 Vmath::Scatr(map.size(), sign.data(), local.data(), map.data(),
820 global.data());
821}
822
824 Array<OneD, NekDouble> &loc) const
825{
827
828 Array<OneD, const int> map = m_patchMapFromPrevLevel->GetNewLevelMap();
830
831 if (global.data() == loc.data())
832 {
833 glo = Array<OneD, NekDouble>(global.size(), global.data());
834 }
835 else
836 {
837 glo = global; // create reference
838 }
839
840 Vmath::Gathr(map.size(), sign.data(), glo.data(), map.data(), loc.data());
841}
842
844 Array<OneD, NekDouble> &global) const
845{
847 Array<OneD, const int> map = m_patchMapFromPrevLevel->GetNewLevelMap();
849
850 if (global.data() == loc.data())
851 {
852 local = Array<OneD, NekDouble>(map.size(), loc.data());
853 }
854 else
855 {
856 local = loc; // create reference
857 }
858
859 // since we are calling mapping from level down from array
860 // the m_numLocaBndCoeffs represents the size of the
861 // boundary elements we need to assemble into
862 Vmath::Zero(m_numLocalCoeffs, global.data(), 1);
863
864 Vmath::Assmb(map.size(), sign.data(), local.data(), map.data(),
865 global.data());
866}
867
872
877
882
887
889{
890 return v_GetNumDirEdges();
891}
892
894{
895 return v_GetNumDirFaces();
896}
897
899{
900 return v_GetNumNonDirEdges();
901}
902
904{
905 return v_GetNumNonDirFaces();
906}
907
912
913std::shared_ptr<AssemblyMap> AssemblyMap::LinearSpaceMap(
914 const ExpList &locexp, GlobalSysSolnType solnType)
915{
916 return v_LinearSpaceMap(locexp, solnType);
917}
918
920{
921 return m_localToGlobalBndMap[i];
922}
923
928
930{
931 return m_signChange;
932}
933
938
943
948
950{
951 if (m_signChange)
952 {
953 return m_localToGlobalBndSign[i];
954 }
955 else
956 {
957 return 1.0;
958 }
959}
960
965
970
975
977{
978 ASSERTL1(i < m_bndCondIDToGlobalTraceID.size(), "Index out of range.");
980}
981
986
988{
989 ASSERTL1(i < m_perbndCondIDToGlobalTraceID.size(), "Index out of range.");
991}
992
997
1002
1007
1009{
1010 return m_numLocalBndCoeffs;
1011}
1012
1017
1019{
1020 return m_numLocalCoeffs;
1021}
1022
1024{
1025 return m_numGlobalCoeffs;
1026}
1027
1029{
1030 return m_systemSingular;
1031}
1032
1034 NekVector<NekDouble> &loc, int offset) const
1035{
1036 GlobalToLocalBnd(global.GetPtr(), loc.GetPtr(), offset);
1037}
1038
1040 NekVector<NekDouble> &loc) const
1041{
1042 GlobalToLocalBnd(global.GetPtr(), loc.GetPtr());
1043}
1044
1047 int offset) const
1048{
1049 ASSERTL1(loc.size() >= m_numLocalBndCoeffs,
1050 "Local vector is not of correct dimension");
1051 ASSERTL1(global.size() >= m_numGlobalBndCoeffs - offset,
1052 "Global vector is not of correct dimension");
1053
1054 // offset input data by length "offset" for Dirichlet boundary conditions.
1056 Vmath::Vcopy(m_numGlobalBndCoeffs - offset, global.data(), 1,
1057 tmp.data() + offset, 1);
1058
1059 if (m_signChange)
1060 {
1062 tmp.data(), m_localToGlobalBndMap.data(), loc.data());
1063 }
1064 else
1065 {
1067 m_localToGlobalBndMap.data(), loc.data());
1068 }
1069}
1070
1072 Array<OneD, NekDouble> &loc) const
1073{
1074 ASSERTL1(loc.size() >= m_numLocalBndCoeffs,
1075 "Local vector is not of correct dimension");
1076 ASSERTL1(global.size() >= m_numGlobalBndCoeffs,
1077 "Global vector is not of correct dimension");
1078
1080 if (global.data() == loc.data())
1081 {
1082 glo = Array<OneD, NekDouble>(m_numLocalBndCoeffs, global.data());
1083 }
1084 else
1085 {
1086 glo = global; // create reference
1087 }
1088
1089 if (m_signChange)
1090 {
1092 glo.data(), m_localToGlobalBndMap.data(), loc.data());
1093 }
1094 else
1095 {
1097 m_localToGlobalBndMap.data(), loc.data());
1098 }
1099}
1100
1102 Array<OneD, NekDouble> &global, int offset,
1103 bool UseComm) const
1104{
1105 ASSERTL1(loc.size() >= m_numLocalBndCoeffs,
1106 "Local vector is not of correct dimension");
1107 ASSERTL1(global.size() >= m_numGlobalBndCoeffs - offset,
1108 "Global vector is not of correct dimension");
1109
1110 // offset input data by length "offset" for Dirichlet boundary conditions.
1112
1113 if (m_signChange)
1114 {
1116 loc.data(), m_localToGlobalBndMap.data(), tmp.data());
1117 }
1118 else
1119 {
1121 m_localToGlobalBndMap.data(), tmp.data());
1122 }
1123
1124 // Ensure each processor has unique value with a max gather.
1125 if (UseComm)
1126 {
1128 }
1129 Vmath::Vcopy(m_numGlobalBndCoeffs - offset, tmp.data() + offset, 1,
1130 global.data(), 1);
1131}
1132
1134 Array<OneD, NekDouble> &global,
1135 bool UseComm) const
1136{
1137 ASSERTL1(loc.size() >= m_numLocalBndCoeffs,
1138 "Local vector is not of correct dimension");
1139 ASSERTL1(global.size() >= m_numGlobalBndCoeffs,
1140 "Global vector is not of correct dimension");
1141
1142 if (m_signChange)
1143 {
1145 loc.data(), m_localToGlobalBndMap.data(), global.data());
1146 }
1147 else
1148 {
1150 m_localToGlobalBndMap.data(), global.data());
1151 }
1152 if (UseComm)
1153 {
1154 Gs::Gather(global, Gs::gs_max, m_bndGsh);
1155 }
1156}
1157
1159 Array<OneD, NekDouble> &locbnd) const
1160{
1161 ASSERTL1(locbnd.size() >= m_numLocalBndCoeffs,
1162 "LocBnd vector is not of correct dimension");
1163 ASSERTL1(local.size() >= m_numLocalCoeffs,
1164 "Local vector is not of correct dimension");
1165
1167 locbnd.data());
1168}
1169
1171 Array<OneD, NekDouble> &locint) const
1172{
1174 "Locint vector is not of correct dimension");
1175 ASSERTL1(local.size() >= m_numLocalCoeffs,
1176 "Local vector is not of correct dimension");
1177
1179 m_localToLocalIntMap.data(), locint.data());
1180}
1181
1183 Array<OneD, NekDouble> &local) const
1184{
1185 ASSERTL1(locbnd.size() >= m_numLocalBndCoeffs,
1186 "LocBnd vector is not of correct dimension");
1187 ASSERTL1(local.size() >= m_numLocalCoeffs,
1188 "Local vector is not of correct dimension");
1189
1190 Vmath::Scatr(m_numLocalBndCoeffs, locbnd.data(),
1191 m_localToLocalBndMap.data(), local.data());
1192}
1193
1195 Array<OneD, NekDouble> &local) const
1196{
1198 "LocBnd vector is not of correct dimension");
1199 ASSERTL1(local.size() >= m_numLocalCoeffs,
1200 "Local vector is not of correct dimension");
1201
1203 m_localToLocalIntMap.data(), local.data());
1204}
1205
1207 NekVector<NekDouble> &global, int offset) const
1208{
1209 AssembleBnd(loc.GetPtr(), global.GetPtr(), offset);
1210}
1211
1213 NekVector<NekDouble> &global) const
1214{
1215 AssembleBnd(loc.GetPtr(), global.GetPtr());
1216}
1217
1219 Array<OneD, NekDouble> &global, int offset) const
1220{
1221 ASSERTL1(loc.size() >= m_numLocalBndCoeffs,
1222 "Local array is not of correct dimension");
1223 ASSERTL1(global.size() >= m_numGlobalBndCoeffs - offset,
1224 "Global array is not of correct dimension");
1226
1227 if (m_signChange)
1228 {
1230 loc.data(), m_localToGlobalBndMap.data(), tmp.data());
1231 }
1232 else
1233 {
1235 m_localToGlobalBndMap.data(), tmp.data());
1236 }
1238 Vmath::Vcopy(m_numGlobalBndCoeffs - offset, tmp.data() + offset, 1,
1239 global.data(), 1);
1240}
1241
1243 Array<OneD, NekDouble> &global) const
1244{
1245 ASSERTL1(loc.size() >= m_numLocalBndCoeffs,
1246 "Local vector is not of correct dimension");
1247 ASSERTL1(global.size() >= m_numGlobalBndCoeffs,
1248 "Global vector is not of correct dimension");
1249
1250 Vmath::Zero(m_numGlobalBndCoeffs, global.data(), 1);
1251
1252 if (m_signChange)
1253 {
1255 loc.data(), m_localToGlobalBndMap.data(), global.data());
1256 }
1257 else
1258 {
1260 m_localToGlobalBndMap.data(), global.data());
1261 }
1262 UniversalAssembleBnd(global);
1263}
1264
1266{
1267 ASSERTL1(pGlobal.size() >= m_numGlobalBndCoeffs, "Wrong size.");
1268 Gs::Gather(pGlobal, Gs::gs_add, m_bndGsh);
1269}
1270
1275
1277 int offset) const
1278{
1279 Array<OneD, NekDouble> tmp(offset);
1280 if (offset > 0)
1281 {
1282 Vmath::Vcopy(offset, pGlobal, 1, tmp, 1);
1283 }
1284 UniversalAssembleBnd(pGlobal);
1285 if (offset > 0)
1286 {
1287 Vmath::Vcopy(offset, tmp, 1, pGlobal, 1);
1288 }
1289}
1290
1295
1300
1302{
1303 return m_staticCondLevel;
1304}
1305
1307{
1308 return m_numPatches;
1309}
1310
1316
1322
1327
1332
1334{
1335 return !(m_nextLevelLocalToGlobalMap.get());
1336}
1337
1342
1344{
1345 return m_preconType;
1346}
1347
1349{
1350 return m_isAbsoluteTolerance;
1351}
1352
1354{
1355 return m_successiveRHS;
1356}
1357
1359{
1360 return m_linSysIterSolver;
1361}
1362
1365{
1366 ASSERTL1(loc.size() >= m_numLocalBndCoeffs,
1367 "Local vector is not of correct dimension");
1368 ASSERTL1(global.size() >= m_numGlobalBndCoeffs,
1369 "Global vector is not of correct dimension");
1370
1371 Vmath::Gathr(m_numLocalBndCoeffs, global.data(),
1372 m_localToGlobalBndMap.data(), loc.data());
1373}
1374
1375void AssemblyMap::PrintStats(std::ostream &out, std::string variable,
1376 bool printHeader) const
1377{
1378 LibUtilities::CommSharedPtr vRowComm = m_session->GetComm()->GetRowComm();
1379 bool isRoot = m_session->GetComm()->IsParallelInTime()
1380 ? m_session->GetComm()->GetRank() == 0
1381 : vRowComm->GetRank() == 0;
1382 int n = vRowComm->GetSize();
1383 int i;
1384
1385 // Determine number of global degrees of freedom.
1386 int globBndCnt = 0, globDirCnt = 0;
1387
1388 for (i = 0; i < m_numGlobalBndCoeffs; ++i)
1389 {
1391 {
1392 globBndCnt++;
1393
1395 {
1396 globDirCnt++;
1397 }
1398 }
1399 }
1400
1401 int globCnt = m_numGlobalCoeffs - m_numGlobalBndCoeffs + globBndCnt;
1402
1403 // Calculate maximum valency
1406
1407 Vmath::Assmb(m_numLocalBndCoeffs, tmpLoc.data(),
1408 m_localToGlobalBndMap.data(), tmpGlob.data());
1409 UniversalAssembleBnd(tmpGlob);
1410
1411 int totGlobDof = globCnt;
1412 int totGlobBndDof = globBndCnt;
1413 int totGlobDirDof = globDirCnt;
1414 int totLocalDof = m_numLocalCoeffs;
1415 int totLocalBndDof = m_numLocalBndCoeffs;
1416 int totLocalDirDof = m_numLocalDirBndCoeffs;
1417
1418 int meanValence = 0;
1419 int maxValence = 0;
1420 int minValence = 10000000;
1421 for (int i = 0; i < m_numGlobalBndCoeffs; ++i)
1422 {
1424 {
1425 continue;
1426 }
1427
1428 if (tmpGlob[i] > maxValence)
1429 {
1430 maxValence = tmpGlob[i];
1431 }
1432 if (tmpGlob[i] < minValence)
1433 {
1434 minValence = tmpGlob[i];
1435 }
1436 meanValence += tmpGlob[i];
1437 }
1438
1439 vRowComm->AllReduce(maxValence, LibUtilities::ReduceMax);
1440 vRowComm->AllReduce(minValence, LibUtilities::ReduceMin);
1441 vRowComm->AllReduce(meanValence, LibUtilities::ReduceSum);
1442 vRowComm->AllReduce(totGlobDof, LibUtilities::ReduceSum);
1443 vRowComm->AllReduce(totGlobBndDof, LibUtilities::ReduceSum);
1444 vRowComm->AllReduce(totGlobDirDof, LibUtilities::ReduceSum);
1445 vRowComm->AllReduce(totLocalDof, LibUtilities::ReduceSum);
1446 vRowComm->AllReduce(totLocalBndDof, LibUtilities::ReduceSum);
1447 vRowComm->AllReduce(totLocalDirDof, LibUtilities::ReduceSum);
1448
1449 meanValence /= totGlobBndDof;
1450
1451 if (isRoot)
1452 {
1453 if (printHeader)
1454 {
1455 out << "Assembly map statistics for field " << variable << ":"
1456 << endl;
1457 }
1458
1459 out << " - Number of local/global dof : " << totLocalDof
1460 << " " << totGlobDof << endl;
1461 out << " - Number of local/global boundary dof : " << totLocalBndDof
1462 << " " << totGlobBndDof << endl;
1463 out << " - Number of local/global Dirichlet dof : " << totLocalDirDof
1464 << " " << totGlobDirDof << endl;
1465 out << " - dof valency (min/max/mean) : " << minValence
1466 << " " << maxValence << " " << meanValence << endl;
1467
1468 if (n > 1)
1469 {
1470 NekDouble mean = m_numLocalCoeffs, mean2 = mean * mean;
1471 NekDouble minval = mean, maxval = mean;
1473
1474 for (i = 1; i < n; ++i)
1475 {
1476 vRowComm->Recv(i, tmp);
1477 mean += tmp[0];
1478 mean2 += tmp[0] * tmp[0];
1479
1480 if (tmp[0] > maxval)
1481 {
1482 maxval = tmp[0];
1483 }
1484 if (tmp[0] < minval)
1485 {
1486 minval = tmp[0];
1487 }
1488 }
1489
1490 if (maxval > 0.1)
1491 {
1492 out << " - Local dof dist. (min/max/mean/dev) : " << minval
1493 << " " << maxval << " " << (mean / n) << " "
1494 << sqrt(mean2 / n - mean * mean / n / n) << endl;
1495 }
1496
1497 vRowComm->Block();
1498
1499 mean = minval = maxval = m_numLocalBndCoeffs;
1500 mean2 = mean * mean;
1501
1502 for (i = 1; i < n; ++i)
1503 {
1504 vRowComm->Recv(i, tmp);
1505 mean += tmp[0];
1506 mean2 += tmp[0] * tmp[0];
1507
1508 if (tmp[0] > maxval)
1509 {
1510 maxval = tmp[0];
1511 }
1512 if (tmp[0] < minval)
1513 {
1514 minval = tmp[0];
1515 }
1516 }
1517
1518 out << " - Local bnd dof dist. (min/max/mean/dev) : " << minval
1519 << " " << maxval << " " << (mean / n) << " "
1520 << sqrt(mean2 / n - mean * mean / n / n) << endl;
1521 }
1522 }
1523 else
1524 {
1526 tmp[0] = m_numLocalCoeffs;
1527 vRowComm->Send(0, tmp);
1528 vRowComm->Block();
1529 tmp[0] = m_numLocalBndCoeffs;
1530 vRowComm->Send(0, tmp);
1531 }
1532
1533 // Either we have no more levels in the static condensation, or we
1534 // are not multi-level.
1536 {
1537 return;
1538 }
1539
1540 int level = 2;
1542 while (tmp->m_nextLevelLocalToGlobalMap)
1543 {
1544 tmp = tmp->m_nextLevelLocalToGlobalMap;
1545 ++level;
1546 }
1547
1548 // Print out multi-level static condensation information.
1549 if (n > 1)
1550 {
1551 if (isRoot)
1552 {
1553 NekDouble mean = level, mean2 = mean * mean;
1554 int minval = level, maxval = level;
1555
1556 Array<OneD, NekDouble> tmpRecv(1);
1557 for (i = 1; i < n; ++i)
1558 {
1559 vRowComm->Recv(i, tmpRecv);
1560 mean += tmpRecv[0];
1561 mean2 += tmpRecv[0] * tmpRecv[0];
1562
1563 if (tmpRecv[0] > maxval)
1564 {
1565 maxval = (int)(tmpRecv[0] + 0.5);
1566 }
1567 if (tmpRecv[0] < minval)
1568 {
1569 minval = (int)(tmpRecv[0] + 0.5);
1570 }
1571 }
1572
1573 out << " - M-level sc. dist. (min/max/mean/dev) : " << minval
1574 << " " << maxval << " " << (mean / n) << " "
1575 << sqrt(mean2 / n - mean * mean / n / n) << endl;
1576 }
1577 else
1578 {
1579 Array<OneD, NekDouble> tmpSend(1);
1580 tmpSend[0] = level;
1581 vRowComm->Send(0, tmpSend);
1582 }
1583 }
1584 else
1585 {
1586 if (isRoot)
1587 {
1588 out << " - Number of static cond. levels : " << level
1589 << endl;
1590 }
1591 }
1592
1593 if (isRoot)
1594 {
1595 out << "Stats at lowest static cond. level:" << endl;
1596 }
1597 tmp->PrintStats(out, variable, false);
1598}
1599} // namespace Nektar::MultiRegions
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
#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.
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.
Array< OneD, int > m_bndCondCoeffsToLocalTraceMap
Integer map of bnd cond coeff to local trace coeff.
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.
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.
Array< OneD, int > m_perbndCondIDToGlobalTraceID
Integer map of rotational periodic bnd cond trace number to global trace number.
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.
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.
Array< OneD, int > m_localToLocalIntMap
Integer map of local boundary coeffs to local interior system numbering.
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMapUnique()
int m_numGlobalCoeffs
Total number of global coefficients.
std::string m_linSysIterSolver
Iterative solver: Conjugate Gradient, GMRES.
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.
void PatchAssemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
int m_successiveRHS
sucessive RHS for iterative solver
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.
void UniversalAbsMaxBnd(Array< OneD, NekDouble > &bndvals)
Array< OneD, unsigned int > m_numLocalBndCoeffsPerPatch
The number of bnd dofs per patch.
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.
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.
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.
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.
int m_bndSystemBandWidth
The bandwith of the global bnd system.
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.
bool GetSingularSystem() const
Retrieves if the system is singular (true) or not (false)
std::string m_variable
Variable string identifier.
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.
LibUtilities::CommSharedPtr GetComm()
Retrieves the communicator.
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
The number of int dofs per patch.
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.
Array< OneD, int > m_localToGlobalBndMap
Integer map of local coeffs to global Boundary Dofs.
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.
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.
Array< OneD, int > m_globalToUniversalBndMapUnique
Integer map of unique process coeffs to universal space (signed)
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.
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.
Array< OneD, NekDouble > m_bndCondCoeffsToLocalCoeffsSign
Integer map of sign of bnd cond coeffs to local coefficients.
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.
const Array< OneD, const int > & GetPerBndCondIDToGlobalTraceID()
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.
int GetNumLocalBndCoeffs() const
Returns the total number of local boundary coefficients.
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Array< OneD, int > m_bndCondIDToGlobalTraceID
Integer map of bnd cond trace number to global trace number.
Base class for all multi-elemental spectral/hp expansions.
Definition ExpList.h:99
Array< OneD, DataType > & GetPtr()
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.
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:290