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#include <boost/core/ignore_unused.hpp>
37
39
40using namespace std;
41
42namespace Nektar
43{
44namespace MultiRegions
45{
46/**
47 * @class AssemblyMap
48 * This class acts as a base class for constructing mappings between
49 * local, global and boundary degrees of freedom. It holds the storage
50 * for the maps and provides the accessors needed to retrieve them.
51 *
52 * There are two derived classes: AssemblyMapCG and
53 * AssemblyMapDG. These perform the actual construction of the
54 * maps within their specific contexts.
55 *
56 */
57
58/// Rounds a double precision number to an integer.
60{
61 return int(x > 0.0 ? x + 0.5 : x - 0.5);
62}
63
64/// Rounds an array of double precision numbers to integers.
66 Array<OneD, int> outarray)
67{
68 int size = inarray.size();
69 ASSERTL1(outarray.size() >= size, "Array sizes not compatible");
70
71 NekDouble x;
72 for (int i = 0; i < size; i++)
73 {
74 x = inarray[i];
75 outarray[i] = int(x > 0.0 ? x + 0.5 : x - 0.5);
76 }
77}
78
79/**
80 * Initialises an empty mapping.
81 */
83 : m_session(), m_comm(), m_hash(0), m_numLocalBndCoeffs(0),
84 m_numGlobalBndCoeffs(0), m_numLocalDirBndCoeffs(0),
85 m_numGlobalDirBndCoeffs(0), m_solnType(eNoSolnType),
86 m_bndSystemBandWidth(0), m_successiveRHS(0),
87 m_linSysIterSolver("ConjugateGradient"), m_gsh(0), m_bndGsh(0)
88{
89}
90
93 const std::string variable)
94 : m_session(pSession), m_comm(comm), m_variable(variable), m_hash(0),
95 m_numLocalBndCoeffs(0), m_numGlobalBndCoeffs(0),
96 m_numLocalDirBndCoeffs(0), m_numGlobalDirBndCoeffs(0),
97 m_bndSystemBandWidth(0), m_successiveRHS(0),
98 m_linSysIterSolver("ConjugateGradient"), m_gsh(0), m_bndGsh(0)
99{
100 // Default value from Solver Info
101 m_solnType =
102 pSession->GetSolverInfoAsEnum<GlobalSysSolnType>("GlobalSysSoln");
103
104 // Override values with data from GlobalSysSolnInfo section
105 if (pSession->DefinesGlobalSysSolnInfo(variable, "GlobalSysSoln"))
106 {
107 std::string sysSoln =
108 pSession->GetGlobalSysSolnInfo(variable, "GlobalSysSoln");
109 m_solnType = pSession->GetValueAsEnum<GlobalSysSolnType>(
110 "GlobalSysSoln", sysSoln);
111 }
112
113 // Set up preconditioner with SysSoln as override then solverinfo otherwise
114 // set a default as diagonal
115 if (pSession->DefinesGlobalSysSolnInfo(variable, "Preconditioner"))
116 {
118 pSession->GetGlobalSysSolnInfo(variable, "Preconditioner");
119 }
120 else if (pSession->DefinesSolverInfo("Preconditioner"))
121 {
122 m_preconType = pSession->GetSolverInfo("Preconditioner");
123 }
124 else
125 { // Possibly preconditioner is default registered in
126 // GlobLinSysIterative.cpp
127 m_preconType = "Diagonal";
128 }
129
130 if (pSession->DefinesGlobalSysSolnInfo(variable,
131 "IterativeSolverTolerance"))
132 {
133 m_iterativeTolerance = boost::lexical_cast<NekDouble>(
134 pSession->GetGlobalSysSolnInfo(variable, "IterativeSolverTolerance")
135 .c_str());
136 }
137 else
138 {
139 pSession->LoadParameter("IterativeSolverTolerance",
142 }
143
144 if (pSession->DefinesGlobalSysSolnInfo(variable, "AbsoluteTolerance"))
145 {
146 std::string abstol =
147 pSession->GetGlobalSysSolnInfo(variable, "AbsoluteTolerance");
149 boost::iequals(boost::to_upper_copy(abstol), "TRUE");
150 }
151 else
152 {
153 m_isAbsoluteTolerance = false;
154 }
155
156 if (pSession->DefinesGlobalSysSolnInfo(variable, "MaxIterations"))
157 {
158 m_maxIterations = boost::lexical_cast<int>(
159 pSession->GetGlobalSysSolnInfo(variable, "MaxIterations").c_str());
160 }
161 else
162 {
163 pSession->LoadParameter("MaxIterations", m_maxIterations, 5000);
164 }
165
166 if (pSession->DefinesGlobalSysSolnInfo(variable, "SuccessiveRHS"))
167 {
168 m_successiveRHS = boost::lexical_cast<int>(
169 pSession->GetGlobalSysSolnInfo(variable, "SuccessiveRHS").c_str());
170 }
171 else
172 {
173 pSession->LoadParameter("SuccessiveRHS", m_successiveRHS, 0);
174 }
175
176 if (pSession->DefinesGlobalSysSolnInfo(variable, "LinSysIterSolver"))
177 {
179 pSession->GetGlobalSysSolnInfo(variable, "LinSysIterSolver");
180 }
181 else if (pSession->DefinesSolverInfo("LinSysIterSolver"))
182 {
183 m_linSysIterSolver = pSession->GetSolverInfo("LinSysIterSolver");
184 }
185 else
186 {
187 m_linSysIterSolver = "ConjugateGradient";
188 }
189}
190
191/**
192 * Create a new level of mapping using the information in
193 * multiLevelGraph and performing the following steps:
194 */
196 AssemblyMap *oldLevelMap,
197 const BottomUpSubStructuredGraphSharedPtr &multiLevelGraph)
198 : m_session(oldLevelMap->m_session), m_comm(oldLevelMap->GetComm()),
199 m_hash(0), m_solnType(oldLevelMap->m_solnType),
200 m_preconType(oldLevelMap->m_preconType),
201 m_maxIterations(oldLevelMap->m_maxIterations),
202 m_iterativeTolerance(oldLevelMap->m_iterativeTolerance),
203 m_successiveRHS(oldLevelMap->m_successiveRHS),
204 m_linSysIterSolver(oldLevelMap->m_linSysIterSolver),
205 m_gsh(oldLevelMap->m_gsh), m_bndGsh(oldLevelMap->m_bndGsh),
206 m_lowestStaticCondLevel(oldLevelMap->m_lowestStaticCondLevel)
207{
208 int i;
209 int j;
210 int cnt;
211
212 //--------------------------------------------------------------
213 // -- Extract information from the input argument
214 int numGlobalBndCoeffsOld = oldLevelMap->GetNumGlobalBndCoeffs();
215 int numGlobalDirBndCoeffsOld = oldLevelMap->GetNumGlobalDirBndCoeffs();
216 int numLocalBndCoeffsOld = oldLevelMap->GetNumLocalBndCoeffs();
217 int numLocalDirBndCoeffsOld = oldLevelMap->GetNumLocalDirBndCoeffs();
218 bool signChangeOld = oldLevelMap->GetSignChange();
219
220 int staticCondLevelOld = oldLevelMap->GetStaticCondLevel();
221 int numPatchesOld = oldLevelMap->GetNumPatches();
222 GlobalSysSolnType solnTypeOld = oldLevelMap->GetGlobalSysSolnType();
223 const Array<OneD, const unsigned int> &numLocalBndCoeffsPerPatchOld =
224 oldLevelMap->GetNumLocalBndCoeffsPerPatch();
225 //--------------------------------------------------------------
226
227 //--------------------------------------------------------------
228 int newLevel = staticCondLevelOld + 1;
229 /** - STEP 1: setup a mask array to determine to which patch
230 * of the new level every patch of the current
231 * level belongs. To do so we make four arrays,
232 * #gloPatchMask, #globHomPatchMask,
233 * #locPatchMask_NekDouble and #locPatchMask.
234 * These arrays are then used to check which local
235 * dofs of the old level belong to which patch of
236 * the new level
237 */
238 Array<OneD, NekDouble> globPatchMask(numGlobalBndCoeffsOld, -1.0);
239 Array<OneD, NekDouble> globHomPatchMask(globPatchMask +
240 numGlobalDirBndCoeffsOld);
241 Array<OneD, NekDouble> locPatchMask_NekDouble(numLocalBndCoeffsOld, -3.0);
242 Array<OneD, int> locPatchMask(numLocalBndCoeffsOld);
243
244 // Fill the array globPatchMask as follows:
245 // - The first part (i.e. the glob bnd dofs) is filled with the
246 // value -1
247 // - The second part (i.e. the glob interior dofs) is numbered
248 // according to the patch it belongs to (i.e. dofs in first block
249 // all are numbered 0, the second block numbered are 1, etc...)
250 multiLevelGraph->MaskPatches(newLevel, globHomPatchMask);
251
252 // Map from Global Dofs to Local Dofs
253 // As a result, we know for each local dof whether
254 // it is mapped to the boundary of the next level, or to which
255 // patch. Based upon this, we can than later associate every patch
256 // of the current level with a patch in the next level.
257 oldLevelMap->GlobalToLocalBndWithoutSign(globPatchMask,
258 locPatchMask_NekDouble);
259
260 // Convert the result to an array of integers rather than doubles
261 RoundNekDoubleToInt(locPatchMask_NekDouble, locPatchMask);
262
263 /** - STEP 2: We calculate how many local bnd dofs of the
264 * old level belong to the boundaries of each patch at
265 * the new level. We need this information to set up the
266 * mapping between different levels.
267 */
268
269 // Retrieve the number of patches at the next level
270 int numPatchesWithIntNew =
271 multiLevelGraph->GetNpatchesWithInterior(newLevel);
272 int numPatchesNew = numPatchesWithIntNew;
273
274 // Allocate memory to store the number of local dofs associated to
275 // each of elemental boundaries of these patches
276 std::map<int, int> numLocalBndCoeffsPerPatchNew;
277 for (int i = 0; i < numPatchesNew; i++)
278 {
279 numLocalBndCoeffsPerPatchNew[i] = 0;
280 }
281
282 int minval;
283 int maxval;
284 int curPatch;
285 for (i = cnt = 0; i < numPatchesOld; i++)
286 {
287 // For every patch at the current level, the mask array
288 // locPatchMask should be filled with
289 // - the same (positive) number for each entry (which will
290 // correspond to the patch at the next level it belongs to)
291 // - the same (positive) number for each entry, except some
292 // entries that are -1 (the enties correspond to -1, will be
293 // mapped to the local boundary of the next level patch given
294 // by the positive number)
295 // - -1 for all entries. In this case, we will make an
296 // additional patch only consisting of boundaries at the next
297 // level
298 minval =
299 *min_element(&locPatchMask[cnt],
300 &locPatchMask[cnt] + numLocalBndCoeffsPerPatchOld[i]);
301 maxval =
302 *max_element(&locPatchMask[cnt],
303 &locPatchMask[cnt] + numLocalBndCoeffsPerPatchOld[i]);
304 ASSERTL0((minval == maxval) || (minval == -1),
305 "These values should never be the same");
306
307 if (maxval == -1)
308 {
309 curPatch = numPatchesNew;
310 numLocalBndCoeffsPerPatchNew[curPatch] = 0;
311 numPatchesNew++;
312 }
313 else
314 {
315 curPatch = maxval;
316 }
317
318 for (j = 0; j < numLocalBndCoeffsPerPatchOld[i]; j++)
319 {
320 ASSERTL0((locPatchMask[cnt] == maxval) ||
321 (locPatchMask[cnt] == minval),
322 "These values should never be the same");
323 if (locPatchMask[cnt] == -1)
324 {
325 ++numLocalBndCoeffsPerPatchNew[curPatch];
326 }
327 cnt++;
328 }
329 }
330
331 // Count how many local dofs of the old level are mapped
332 // to the local boundary dofs of the new level
335 m_numPatches = numLocalBndCoeffsPerPatchNew.size();
338 multiLevelGraph->GetNintDofsPerPatch(newLevel, m_numLocalIntCoeffsPerPatch);
339
340 for (int i = 0; i < m_numPatches; i++)
341 {
343 (unsigned int)numLocalBndCoeffsPerPatchNew[i];
347 }
348
349 // Also initialise some more data members
350 m_solnType = solnTypeOld;
355 "This method should only be called for in "
356 "case of multi-level static condensation.");
357 m_staticCondLevel = newLevel;
358 m_signChange = signChangeOld;
359 m_numLocalDirBndCoeffs = numLocalDirBndCoeffsOld;
360 m_numGlobalDirBndCoeffs = numGlobalDirBndCoeffsOld;
362 multiLevelGraph->GetInteriorOffset(newLevel) + m_numGlobalDirBndCoeffs;
364 multiLevelGraph->GetNumGlobalDofs(newLevel) + m_numGlobalDirBndCoeffs;
366
370
371 // local to bnd map is just a copy
372 for (int i = 0; i < m_numLocalBndCoeffs; ++i)
373 {
375 }
376
377 // local to int map is just a copy plus offset
378 for (int i = m_numLocalBndCoeffs; i < m_numLocalCoeffs; ++i)
379 {
381 }
382
383 if (m_signChange)
384 {
386 }
387
389 MemoryManager<PatchMap>::AllocateSharedPtr(numLocalBndCoeffsOld);
390
395
396 // Set up an offset array that denotes the offset of the local
397 // boundary degrees of freedom of the next level
398 Array<OneD, int> numLocalBndCoeffsPerPatchOffset(m_numPatches + 1, 0);
399 for (int i = 1; i < m_numPatches + 1; i++)
400 {
401 numLocalBndCoeffsPerPatchOffset[i] +=
402 numLocalBndCoeffsPerPatchOffset[i - 1] +
403 numLocalBndCoeffsPerPatchNew[i - 1];
404 }
405
406 int additionalPatchCnt = numPatchesWithIntNew;
407 int newid;
408 int blockid;
409 bool isBndDof;
411 Array<OneD, int> bndDofPerPatchCnt(m_numPatches, 0);
412
413 for (i = cnt = 0; i < numPatchesOld; i++)
414 {
415 minval =
416 *min_element(&locPatchMask[cnt],
417 &locPatchMask[cnt] + numLocalBndCoeffsPerPatchOld[i]);
418 maxval =
419 *max_element(&locPatchMask[cnt],
420 &locPatchMask[cnt] + numLocalBndCoeffsPerPatchOld[i]);
421 ASSERTL0((minval == maxval) || (minval == -1),
422 "These values should never be the same");
423
424 if (maxval == -1)
425 {
426 curPatch = additionalPatchCnt;
427 additionalPatchCnt++;
428 }
429 else
430 {
431 curPatch = maxval;
432 }
433
434 for (j = 0; j < numLocalBndCoeffsPerPatchOld[i]; j++)
435 {
436 ASSERTL0((locPatchMask[cnt] == maxval) ||
437 (locPatchMask[cnt] == minval),
438 "These values should never be the same");
439
440 sign = oldLevelMap->GetLocalToGlobalBndSign(cnt);
441
442 if (locPatchMask[cnt] == -1)
443 {
444 newid = numLocalBndCoeffsPerPatchOffset[curPatch];
445
446 m_localToGlobalBndMap[newid] =
447 oldLevelMap->GetLocalToGlobalBndMap(cnt);
448
449 if (m_signChange)
450 {
452 }
453
454 blockid = bndDofPerPatchCnt[curPatch];
455 isBndDof = true;
456
457 numLocalBndCoeffsPerPatchOffset[curPatch]++;
458 bndDofPerPatchCnt[curPatch]++;
459 }
460 else
461 {
462 newid = oldLevelMap->GetLocalToGlobalBndMap(cnt) -
464
465 blockid =
466 oldLevelMap->GetLocalToGlobalBndMap(cnt) -
468 multiLevelGraph->GetInteriorOffset(newLevel, curPatch);
469 isBndDof = false;
470 }
471
472 sign = isBndDof ? 1.0 : sign;
473
474 m_patchMapFromPrevLevel->SetPatchMap(cnt, curPatch, blockid,
475 isBndDof, sign);
476 cnt++;
477 }
478 }
479
480 // set up local to local mapping from previous to new level
483
485
486 // Postprocess the computed information - Update the old
487 // level with the mapping to new level
488 // oldLevelMap->SetLocalBndToNextLevelMap(oldLocalBndToNewLevelMap,oldLocalBndToNewLevelSign);
489
490 // - Construct the next level mapping object
491 if (m_staticCondLevel < (multiLevelGraph->GetNlevels() - 1))
492 {
495 multiLevelGraph);
496 }
497}
498
500{
501}
502
503/**
504 * The bandwidth calculated corresponds to what is referred to as
505 * half-bandwidth. If the elements of the matrix are designated as
506 * a_ij, it corresponds to the maximum value of |i-j| for non-zero
507 * a_ij. As a result, the value also corresponds to the number of
508 * sub- or super-diagonals.
509 *
510 * The bandwith can be calculated elementally as it corresponds to the
511 * maximal elemental bandwith (i.e. the maximal difference in global
512 * DOF index for every element).
513 *
514 * We here calculate the bandwith of the global boundary system (as
515 * used for static condensation).
516 */
518{
519 int i, j;
520 int cnt = 0;
521 int locSize;
522 int maxId;
523 int minId;
524 int bwidth = -1;
525 for (i = 0; i < m_numPatches; ++i)
526 {
527 locSize = m_numLocalBndCoeffsPerPatch[i];
528 maxId = -1;
529 minId = m_numLocalBndCoeffs + 1;
530 for (j = 0; j < locSize; j++)
531 {
533 {
534 if (m_localToGlobalBndMap[cnt + j] > maxId)
535 {
536 maxId = m_localToGlobalBndMap[cnt + j];
537 }
538
539 if (m_localToGlobalBndMap[cnt + j] < minId)
540 {
541 minId = m_localToGlobalBndMap[cnt + j];
542 }
543 }
544 }
545 bwidth = (bwidth > (maxId - minId)) ? bwidth : (maxId - minId);
546
547 cnt += locSize;
548 }
549
550 m_bndSystemBandWidth = bwidth;
551}
552
554{
555 boost::ignore_unused(i);
556 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
557 return 0;
558}
559
561{
562 boost::ignore_unused(i);
563 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
564 return 0;
565}
566
568{
569 boost::ignore_unused(i);
570 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
571 return 0;
572}
573
575{
576 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
577 static Array<OneD, const int> result;
578 return result;
579}
580
582{
583 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
584 static Array<OneD, const int> result;
585 return result;
586}
587
589{
590 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
591 static Array<OneD, const int> result;
592 return result;
593}
594
596{
597 boost::ignore_unused(i);
598 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
599 return 0.0;
600}
601
603{
604 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
605 static Array<OneD, NekDouble> result;
606 return result;
607}
608
611 bool useComm) const
612{
613 boost::ignore_unused(loc, global, useComm);
614 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
615}
616
618 NekVector<NekDouble> &global,
619 bool useComm) const
620{
621 boost::ignore_unused(loc, global, useComm);
622 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
623}
624
627{
628 boost::ignore_unused(loc, global);
629 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
630}
631
634{
635 boost::ignore_unused(loc, global);
636 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
637}
638
640 Array<OneD, NekDouble> &global) const
641{
642 boost::ignore_unused(loc, global);
643 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
644}
645
647 NekVector<NekDouble> &global) const
648{
649 boost::ignore_unused(loc, global);
650 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
651}
652
654{
655 boost::ignore_unused(pGlobal);
656 // Do nothing here since multi-level static condensation uses a
657 // AssemblyMap and thus will call this routine in serial.
658}
659
661{
662 boost::ignore_unused(pGlobal);
663 // Do nothing here since multi-level static condensation uses a
664 // AssemblyMap and thus will call this routine in serial.
665}
666
668 int offset) const
669{
670 boost::ignore_unused(pGlobal, offset);
671 // Do nothing here since multi-level static condensation uses a
672 // AssemblyMap and thus will call this routine in serial.
673}
674
676{
677 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
678 return 0;
679}
680
682{
683 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
684 return 0;
685}
686
688{
689 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
690 return 0;
691}
692
694{
695 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
696 return 0;
697}
698
700{
701 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
702 return 0;
703}
704
706{
707 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
708 return 0;
709}
710
712{
713 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
714 return 0;
715}
716
718{
719 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
720 return 0;
721}
722
724{
725 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
726 static Array<OneD, const int> result;
727 return result;
728}
729
730std::shared_ptr<AssemblyMap> AssemblyMap::v_LinearSpaceMap(
731 const ExpList &locexp, GlobalSysSolnType solnType)
732{
733 boost::ignore_unused(locexp, solnType);
734 NEKERROR(ErrorUtil::efatal, "Not defined for this type of mapping.");
735 static std::shared_ptr<AssemblyMap> result;
736 return result;
737}
738
740{
741 return m_comm;
742}
743
745{
746 return m_variable;
747}
748
750{
751 return m_hash;
752}
753
755{
756 return v_GetLocalToGlobalMap(i);
757}
758
760{
762}
763
765{
767}
768
770{
771 return v_GetLocalToGlobalMap();
772}
773
775{
777}
778
780{
782}
783
785{
786 return v_GetLocalToGlobalSign(i);
787}
788
790{
791 return v_GetLocalToGlobalSign();
792}
793
796 bool useComm) const
797{
798 v_LocalToGlobal(loc, global, useComm);
799}
800
802 NekVector<NekDouble> &global,
803 bool useComm) const
804{
805 v_LocalToGlobal(loc, global, useComm);
806}
807
810{
811 v_GlobalToLocal(global, loc);
812}
813
816{
817 v_GlobalToLocal(global, loc);
818}
819
821 Array<OneD, NekDouble> &global) const
822{
823 v_Assemble(loc, global);
824}
825
827 NekVector<NekDouble> &global) const
828{
829 v_Assemble(loc, global);
830}
831
833{
834 v_UniversalAssemble(pGlobal);
835}
836
838{
839 v_UniversalAssemble(pGlobal);
840}
841
843 int offset) const
844{
845 v_UniversalAssemble(pGlobal, offset);
846}
847
849 Array<OneD, NekDouble> &global) const
850{
852
853 Array<OneD, const int> map = m_patchMapFromPrevLevel->GetNewLevelMap();
855
856 if (global.data() == loc.data())
857 {
858 local = Array<OneD, NekDouble>(map.size(), loc.data());
859 }
860 else
861 {
862 local = loc; // create reference
863 }
864
865 Vmath::Scatr(map.size(), sign.get(), local.get(), map.get(), global.get());
866}
867
870{
872
873 Array<OneD, const int> map = m_patchMapFromPrevLevel->GetNewLevelMap();
875
876 if (global.data() == loc.data())
877 {
878 glo = Array<OneD, NekDouble>(global.size(), global.data());
879 }
880 else
881 {
882 glo = global; // create reference
883 }
884
885 Vmath::Gathr(map.size(), sign.get(), glo.get(), map.get(), loc.get());
886}
887
889 Array<OneD, NekDouble> &global) const
890{
892 Array<OneD, const int> map = m_patchMapFromPrevLevel->GetNewLevelMap();
894
895 if (global.data() == loc.data())
896 {
897 local = Array<OneD, NekDouble>(map.size(), loc.data());
898 }
899 else
900 {
901 local = loc; // create reference
902 }
903
904 // since we are calling mapping from level down from array
905 // the m_numLocaBndCoeffs represents the size of the
906 // boundary elements we need to assemble into
907 Vmath::Zero(m_numLocalCoeffs, global.get(), 1);
908
909 Vmath::Assmb(map.size(), sign.get(), local.get(), map.get(), global.get());
910}
911
913{
915}
916
918{
920}
921
923{
925}
926
928{
930}
931
933{
934 return v_GetNumDirEdges();
935}
936
938{
939 return v_GetNumDirFaces();
940}
941
943{
944 return v_GetNumNonDirEdges();
945}
946
948{
949 return v_GetNumNonDirFaces();
950}
951
953{
954 return v_GetExtraDirEdges();
955}
956
957std::shared_ptr<AssemblyMap> AssemblyMap::LinearSpaceMap(
958 const ExpList &locexp, GlobalSysSolnType solnType)
959{
960 return v_LinearSpaceMap(locexp, solnType);
961}
962
964{
965 return m_localToGlobalBndMap[i];
966}
967
969{
971}
972
974{
975 return m_signChange;
976}
977
979{
981}
982
984{
986}
987
989{
991}
992
994{
995 if (m_signChange)
996 {
997 return m_localToGlobalBndSign[i];
998 }
999 else
1000 {
1001 return 1.0;
1002 }
1003}
1004
1006{
1008}
1009
1011{
1013}
1014
1016{
1018}
1019
1021{
1022 ASSERTL1(i < m_bndCondIDToGlobalTraceID.size(), "Index out of range.");
1024}
1025
1027{
1029}
1030
1032{
1034}
1035
1037{
1039}
1040
1042{
1043 return m_numLocalBndCoeffs;
1044}
1045
1047{
1048 return m_numGlobalBndCoeffs;
1049}
1050
1052{
1053 return m_numLocalCoeffs;
1054}
1055
1057{
1058 return m_numGlobalCoeffs;
1059}
1060
1062{
1063 return m_systemSingular;
1064}
1065
1067 NekVector<NekDouble> &loc, int offset) const
1068{
1069 GlobalToLocalBnd(global.GetPtr(), loc.GetPtr(), offset);
1070}
1071
1074{
1075 GlobalToLocalBnd(global.GetPtr(), loc.GetPtr());
1076}
1077
1080 int offset) const
1081{
1083 "Local vector is not of correct dimension");
1084 ASSERTL1(global.size() >= m_numGlobalBndCoeffs - offset,
1085 "Global vector is not of correct dimension");
1086
1087 // offset input data by length "offset" for Dirichlet boundary conditions.
1089 Vmath::Vcopy(m_numGlobalBndCoeffs - offset, global.get(), 1,
1090 tmp.get() + offset, 1);
1091
1092 if (m_signChange)
1093 {
1095 tmp.get(), m_localToGlobalBndMap.get(), loc.get());
1096 }
1097 else
1098 {
1100 m_localToGlobalBndMap.get(), loc.get());
1101 }
1102}
1103
1106{
1108 "Local vector is not of correct dimension");
1109 ASSERTL1(global.size() >= m_numGlobalBndCoeffs,
1110 "Global vector is not of correct dimension");
1111
1113 if (global.data() == loc.data())
1114 {
1115 glo = Array<OneD, NekDouble>(m_numLocalBndCoeffs, global.data());
1116 }
1117 else
1118 {
1119 glo = global; // create reference
1120 }
1121
1122 if (m_signChange)
1123 {
1125 glo.get(), m_localToGlobalBndMap.get(), loc.get());
1126 }
1127 else
1128 {
1130 m_localToGlobalBndMap.get(), loc.get());
1131 }
1132}
1133
1135 Array<OneD, NekDouble> &global, int offset,
1136 bool UseComm) const
1137{
1139 "Local vector is not of correct dimension");
1140 ASSERTL1(global.size() >= m_numGlobalBndCoeffs - offset,
1141 "Global vector is not of correct dimension");
1142
1143 // offset input data by length "offset" for Dirichlet boundary conditions.
1145
1146 if (m_signChange)
1147 {
1149 loc.get(), m_localToGlobalBndMap.get(), tmp.get());
1150 }
1151 else
1152 {
1154 m_localToGlobalBndMap.get(), tmp.get());
1155 }
1156
1157 // Ensure each processor has unique value with a max gather.
1158 if (UseComm)
1159 {
1161 }
1162 Vmath::Vcopy(m_numGlobalBndCoeffs - offset, tmp.get() + offset, 1,
1163 global.get(), 1);
1164}
1165
1167 Array<OneD, NekDouble> &global,
1168 bool UseComm) const
1169{
1171 "Local vector is not of correct dimension");
1172 ASSERTL1(global.size() >= m_numGlobalBndCoeffs,
1173 "Global vector is not of correct dimension");
1174
1175 if (m_signChange)
1176 {
1178 loc.get(), m_localToGlobalBndMap.get(), global.get());
1179 }
1180 else
1181 {
1183 m_localToGlobalBndMap.get(), global.get());
1184 }
1185 if (UseComm)
1186 {
1187 Gs::Gather(global, Gs::gs_max, m_bndGsh);
1188 }
1189}
1190
1192 Array<OneD, NekDouble> &locbnd) const
1193{
1194 ASSERTL1(locbnd.size() >= m_numLocalBndCoeffs,
1195 "LocBnd vector is not of correct dimension");
1196 ASSERTL1(local.size() >= m_numLocalCoeffs,
1197 "Local vector is not of correct dimension");
1198
1200 locbnd.get());
1201}
1202
1204 Array<OneD, NekDouble> &locint) const
1205{
1207 "Locint vector is not of correct dimension");
1208 ASSERTL1(local.size() >= m_numLocalCoeffs,
1209 "Local vector is not of correct dimension");
1210
1212 m_localToLocalIntMap.get(), locint.get());
1213}
1214
1216 Array<OneD, NekDouble> &local) const
1217{
1218 ASSERTL1(locbnd.size() >= m_numLocalBndCoeffs,
1219 "LocBnd vector is not of correct dimension");
1220 ASSERTL1(local.size() >= m_numLocalCoeffs,
1221 "Local vector is not of correct dimension");
1222
1224 local.get());
1225}
1226
1228 Array<OneD, NekDouble> &local) const
1229{
1231 "LocBnd vector is not of correct dimension");
1232 ASSERTL1(local.size() >= m_numLocalCoeffs,
1233 "Local vector is not of correct dimension");
1234
1236 m_localToLocalIntMap.get(), local.get());
1237}
1238
1240 NekVector<NekDouble> &global, int offset) const
1241{
1242 AssembleBnd(loc.GetPtr(), global.GetPtr(), offset);
1243}
1244
1246 NekVector<NekDouble> &global) const
1247{
1248 AssembleBnd(loc.GetPtr(), global.GetPtr());
1249}
1250
1252 Array<OneD, NekDouble> &global, int offset) const
1253{
1255 "Local array is not of correct dimension");
1256 ASSERTL1(global.size() >= m_numGlobalBndCoeffs - offset,
1257 "Global array is not of correct dimension");
1259
1260 if (m_signChange)
1261 {
1263 loc.get(), m_localToGlobalBndMap.get(), tmp.get());
1264 }
1265 else
1266 {
1268 m_localToGlobalBndMap.get(), tmp.get());
1269 }
1271 Vmath::Vcopy(m_numGlobalBndCoeffs - offset, tmp.get() + offset, 1,
1272 global.get(), 1);
1273}
1274
1276 Array<OneD, NekDouble> &global) const
1277{
1279 "Local vector is not of correct dimension");
1280 ASSERTL1(global.size() >= m_numGlobalBndCoeffs,
1281 "Global vector is not of correct dimension");
1282
1283 Vmath::Zero(m_numGlobalBndCoeffs, global.get(), 1);
1284
1285 if (m_signChange)
1286 {
1288 loc.get(), m_localToGlobalBndMap.get(), global.get());
1289 }
1290 else
1291 {
1293 m_localToGlobalBndMap.get(), global.get());
1294 }
1295 UniversalAssembleBnd(global);
1296}
1297
1299{
1300 ASSERTL1(pGlobal.size() >= m_numGlobalBndCoeffs, "Wrong size.");
1301 Gs::Gather(pGlobal, Gs::gs_add, m_bndGsh);
1302}
1303
1305{
1306 UniversalAssembleBnd(pGlobal.GetPtr());
1307}
1308
1310 int offset) const
1311{
1312 Array<OneD, NekDouble> tmp(offset);
1313 if (offset > 0)
1314 Vmath::Vcopy(offset, pGlobal, 1, tmp, 1);
1315 UniversalAssembleBnd(pGlobal);
1316 if (offset > 0)
1317 Vmath::Vcopy(offset, tmp, 1, pGlobal, 1);
1318}
1319
1321{
1323}
1324
1326{
1327 return m_bndSystemBandWidth;
1328}
1329
1331{
1332 return m_staticCondLevel;
1333}
1334
1336{
1337 return m_numPatches;
1338}
1339
1342{
1344}
1345
1348{
1350}
1351
1353{
1355}
1356
1358{
1360}
1361
1363{
1364 return !(m_nextLevelLocalToGlobalMap.get());
1365}
1366
1368{
1369 return m_solnType;
1370}
1371
1373{
1374 return m_preconType;
1375}
1376
1378{
1379 return m_iterativeTolerance;
1380}
1381
1383{
1384 return m_isAbsoluteTolerance;
1385}
1386
1388{
1389 return m_maxIterations;
1390}
1391
1393{
1394 return m_successiveRHS;
1395}
1396
1398{
1399 return m_linSysIterSolver;
1400}
1401
1404{
1406 "Local vector is not of correct dimension");
1407 ASSERTL1(global.size() >= m_numGlobalBndCoeffs,
1408 "Global vector is not of correct dimension");
1409
1411 loc.get());
1412}
1413
1414void AssemblyMap::PrintStats(std::ostream &out, std::string variable,
1415 bool printHeader) const
1416{
1417 LibUtilities::CommSharedPtr vRowComm = m_session->GetComm()->GetRowComm();
1418 bool isRoot = vRowComm->GetRank() == 0;
1419 int n = vRowComm->GetSize();
1420 int i;
1421
1422 // Determine number of global degrees of freedom.
1423 int globBndCnt = 0, globDirCnt = 0;
1424
1425 for (i = 0; i < m_numGlobalBndCoeffs; ++i)
1426 {
1428 {
1429 globBndCnt++;
1430
1432 {
1433 globDirCnt++;
1434 }
1435 }
1436 }
1437
1438 int globCnt = m_numGlobalCoeffs - m_numGlobalBndCoeffs + globBndCnt;
1439
1440 // Calculate maximum valency
1443
1445 tmpGlob.get());
1446 UniversalAssembleBnd(tmpGlob);
1447
1448 int totGlobDof = globCnt;
1449 int totGlobBndDof = globBndCnt;
1450 int totGlobDirDof = globDirCnt;
1451 int totLocalDof = m_numLocalCoeffs;
1452 int totLocalBndDof = m_numLocalBndCoeffs;
1453 int totLocalDirDof = m_numLocalDirBndCoeffs;
1454
1455 int meanValence = 0;
1456 int maxValence = 0;
1457 int minValence = 10000000;
1458 for (int i = 0; i < m_numGlobalBndCoeffs; ++i)
1459 {
1461 {
1462 continue;
1463 }
1464
1465 if (tmpGlob[i] > maxValence)
1466 {
1467 maxValence = tmpGlob[i];
1468 }
1469 if (tmpGlob[i] < minValence)
1470 {
1471 minValence = tmpGlob[i];
1472 }
1473 meanValence += tmpGlob[i];
1474 }
1475
1476 vRowComm->AllReduce(maxValence, LibUtilities::ReduceMax);
1477 vRowComm->AllReduce(minValence, LibUtilities::ReduceMin);
1478 vRowComm->AllReduce(meanValence, LibUtilities::ReduceSum);
1479 vRowComm->AllReduce(totGlobDof, LibUtilities::ReduceSum);
1480 vRowComm->AllReduce(totGlobBndDof, LibUtilities::ReduceSum);
1481 vRowComm->AllReduce(totGlobDirDof, LibUtilities::ReduceSum);
1482 vRowComm->AllReduce(totLocalDof, LibUtilities::ReduceSum);
1483 vRowComm->AllReduce(totLocalBndDof, LibUtilities::ReduceSum);
1484 vRowComm->AllReduce(totLocalDirDof, LibUtilities::ReduceSum);
1485
1486 meanValence /= totGlobBndDof;
1487
1488 if (isRoot)
1489 {
1490 if (printHeader)
1491 {
1492 out << "Assembly map statistics for field " << variable << ":"
1493 << endl;
1494 }
1495
1496 out << " - Number of local/global dof : " << totLocalDof
1497 << " " << totGlobDof << endl;
1498 out << " - Number of local/global boundary dof : " << totLocalBndDof
1499 << " " << totGlobBndDof << endl;
1500 out << " - Number of local/global Dirichlet dof : " << totLocalDirDof
1501 << " " << totGlobDirDof << endl;
1502 out << " - dof valency (min/max/mean) : " << minValence
1503 << " " << maxValence << " " << meanValence << endl;
1504
1505 if (n > 1)
1506 {
1507 NekDouble mean = m_numLocalCoeffs, mean2 = mean * mean;
1508 NekDouble minval = mean, maxval = mean;
1510
1511 for (i = 1; i < n; ++i)
1512 {
1513 vRowComm->Recv(i, tmp);
1514 mean += tmp[0];
1515 mean2 += tmp[0] * tmp[0];
1516
1517 if (tmp[0] > maxval)
1518 {
1519 maxval = tmp[0];
1520 }
1521 if (tmp[0] < minval)
1522 {
1523 minval = tmp[0];
1524 }
1525 }
1526
1527 if (maxval > 0.1)
1528 {
1529 out << " - Local dof dist. (min/max/mean/dev) : " << minval
1530 << " " << maxval << " " << (mean / n) << " "
1531 << sqrt(mean2 / n - mean * mean / n / n) << endl;
1532 }
1533
1534 vRowComm->Block();
1535
1536 mean = minval = maxval = m_numLocalBndCoeffs;
1537 mean2 = mean * mean;
1538
1539 for (i = 1; i < n; ++i)
1540 {
1541 vRowComm->Recv(i, tmp);
1542 mean += tmp[0];
1543 mean2 += tmp[0] * tmp[0];
1544
1545 if (tmp[0] > maxval)
1546 {
1547 maxval = tmp[0];
1548 }
1549 if (tmp[0] < minval)
1550 {
1551 minval = tmp[0];
1552 }
1553 }
1554
1555 out << " - Local bnd dof dist. (min/max/mean/dev) : " << minval
1556 << " " << maxval << " " << (mean / n) << " "
1557 << sqrt(mean2 / n - mean * mean / n / n) << endl;
1558 }
1559 }
1560 else
1561 {
1563 tmp[0] = m_numLocalCoeffs;
1564 vRowComm->Send(0, tmp);
1565 vRowComm->Block();
1566 tmp[0] = m_numLocalBndCoeffs;
1567 vRowComm->Send(0, tmp);
1568 }
1569
1570 // Either we have no more levels in the static condensation, or we
1571 // are not multi-level.
1573 {
1574 return;
1575 }
1576
1577 int level = 2;
1579 while (tmp->m_nextLevelLocalToGlobalMap)
1580 {
1581 tmp = tmp->m_nextLevelLocalToGlobalMap;
1582 ++level;
1583 }
1584
1585 // Print out multi-level static condensation information.
1586 if (n > 1)
1587 {
1588 if (isRoot)
1589 {
1590 NekDouble mean = level, mean2 = mean * mean;
1591 int minval = level, maxval = level;
1592
1593 Array<OneD, NekDouble> tmpRecv(1);
1594 for (i = 1; i < n; ++i)
1595 {
1596 vRowComm->Recv(i, tmpRecv);
1597 mean += tmpRecv[0];
1598 mean2 += tmpRecv[0] * tmpRecv[0];
1599
1600 if (tmpRecv[0] > maxval)
1601 {
1602 maxval = (int)(tmpRecv[0] + 0.5);
1603 }
1604 if (tmpRecv[0] < minval)
1605 {
1606 minval = (int)(tmpRecv[0] + 0.5);
1607 }
1608 }
1609
1610 out << " - M-level sc. dist. (min/max/mean/dev) : " << minval
1611 << " " << maxval << " " << (mean / n) << " "
1612 << sqrt(mean2 / n - mean * mean / n / n) << endl;
1613 }
1614 else
1615 {
1616 Array<OneD, NekDouble> tmpSend(1);
1617 tmpSend[0] = level;
1618 vRowComm->Send(0, tmpSend);
1619 }
1620 }
1621 else
1622 {
1623 out << " - Number of static cond. levels : " << level << endl;
1624 }
1625
1626 if (isRoot)
1627 {
1628 out << "Stats at lowest static cond. level:" << endl;
1629 }
1630 tmp->PrintStats(out, variable, false);
1631}
1632} // namespace MultiRegions
1633} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:49
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:58
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:405
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:367
Array< OneD, int > m_bndCondCoeffsToLocalTraceMap
Integer map of bnd cond coeff to local trace coeff.
Definition: AssemblyMap.h:396
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:392
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:381
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:520
int m_maxIterations
Maximum iterations for iterative solver.
Definition: AssemblyMap.h:413
Array< OneD, int > m_localToLocalIntMap
Integer map of local boundary coeffs to local interior system numbering.
Definition: AssemblyMap.h:390
virtual const Array< OneD, const int > & v_GetGlobalToUniversalMapUnique()
NekDouble m_iterativeTolerance
Tolerance for iterative solver.
Definition: AssemblyMap.h:416
NekDouble GetIterativeTolerance() const
int m_numGlobalCoeffs
Total number of global coefficients.
Definition: AssemblyMap.h:378
std::string m_linSysIterSolver
Iterative solver: Conjugate Gradient, GMRES.
Definition: AssemblyMap.h:423
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:400
void PatchAssemble(const Array< OneD, const NekDouble > &loc, Array< OneD, NekDouble > &global) const
int m_successiveRHS
sucessive RHS for iterative solver
Definition: AssemblyMap.h:420
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:386
void UniversalAbsMaxBnd(Array< OneD, NekDouble > &bndvals)
Array< OneD, unsigned int > m_numLocalBndCoeffsPerPatch
The number of bnd dofs per patch.
Definition: AssemblyMap.h:436
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:336
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:348
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:443
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:432
int m_bndSystemBandWidth
The bandwith of the global bnd system.
Definition: AssemblyMap.h:407
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:352
bool GetSingularSystem() const
Retrieves if the system is singular (true) or not (false)
std::string m_variable
Variable string identifier.
Definition: AssemblyMap.h:342
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:354
LibUtilities::CommSharedPtr GetComm()
Retrieves the communicator.
Array< OneD, unsigned int > m_numLocalIntCoeffsPerPatch
The number of int dofs per patch.
Definition: AssemblyMap.h:438
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:356
Array< OneD, int > m_localToGlobalBndMap
Integer map of local coeffs to global Boundary Dofs.
Definition: AssemblyMap.h:384
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:428
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:82
Array< OneD, int > m_globalToUniversalBndMapUnique
Integer map of unique process coeffs to universal space (signed)
Definition: AssemblyMap.h:402
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:339
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:388
Array< OneD, NekDouble > m_bndCondCoeffsToLocalCoeffsSign
Integer map of sign of bnd cond coeffs to local coefficients.
Definition: AssemblyMap.h:394
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:434
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:410
int GetNumLocalBndCoeffs() const
Returns the total number of local boundary coefficients.
int m_numGlobalBndCoeffs
Total number of global boundary coefficients.
Definition: AssemblyMap.h:350
Array< OneD, int > m_bndCondIDToGlobalTraceID
Integer map of bnd cond trace number to global trace number.
Definition: AssemblyMap.h:398
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:102
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:270
@ gs_amax
Definition: GsLib.hpp:66
@ gs_add
Definition: GsLib.hpp:62
@ gs_max
Definition: GsLib.hpp:65
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:57
@ eNoSolnType
No Solution type specified.
std::shared_ptr< BottomUpSubStructuredGraph > BottomUpSubStructuredGraphSharedPtr
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:52
std::shared_ptr< PatchMap > PatchMapSharedPtr
int RoundNekDoubleToInt(NekDouble x)
Rounds a double precision number to an integer.
Definition: AssemblyMap.cpp:59
static const NekDouble kNekIterativeTol
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble
void Scatr(int n, const T *x, const int *y, T *z)
Scatter vector z[y[i]] = x[i].
Definition: Vmath.cpp:817
void Gathr(int n, const T *sign, const T *x, const int *y, T *z)
Gather vector z[i] = sign[i]*x[y[i]].
Definition: Vmath.cpp:800
void Assmb(int n, const T *x, const int *y, T *z)
Assemble z[y[i]] += x[i]; z should be zero'd first.
Definition: Vmath.cpp:857
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294