Nektar++
PreconditionerBlock.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: PreconditionerBlock.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: Block Preconditioner definition
32//
33///////////////////////////////////////////////////////////////////////////////
34
36#include <LocalRegions/SegExp.h>
42#include <cmath>
43
44using namespace std;
45
47{
48/**
49 * Registers the class with the Factory.
50 */
53 "Block", PreconditionerBlock::create, "Block Diagonal Preconditioning");
54/**
55 * @class Block Preconditioner
56 *
57 * This class implements Block preconditioning for the conjugate
58 * gradient matrix solver.
59 */
60
62 const std::shared_ptr<GlobalLinSys> &plinsys,
63 const AssemblyMapSharedPtr &pLocToGloMap)
64 : Preconditioner(plinsys, pLocToGloMap)
65{
66}
67
69{
70 GlobalSysSolnType solverType = m_locToGloMap.lock()->GetGlobalSysSolnType();
71 m_isFull =
72 (solverType == eIterativeFull) || (solverType == ePETScFullMatrix);
73
76 "Solver type not valid");
77}
78
80{
81 GlobalLinSysKey key = m_linsys.lock()->GetKey();
82
83 // Different setup for HDG and CG.
85 {
87 }
88 else
89 {
91 }
92}
93
94/**
95 * \brief Construct a block preconditioner from \f$\mathbf{S}_{1}\f$ for
96 * the continuous Galerkin system.
97 *
98 * The preconditioner for the statically-condensed system is defined as:
99 *
100 * \f[\mathbf{M}^{-1}=\left[\begin{array}{ccc}
101 * \mathrm{Diag}[(\mathbf{S_{1}})_{vv}] & & \\
102 * & (\mathbf{S}_{1})_{eb} & \\
103 * & & (\mathbf{S}_{1})_{fb} \end{array}\right] \f]
104 *
105 * where \f$\mathbf{S}_{1}\f$ is the local Schur complement matrix for
106 * each element and the subscript denotes the portion of the Schur
107 * complement associated with the vertex (vv), edge (eb) and face
108 * blocks (fb) respectively.
109 *
110 * In case of the full system \f$\mathbf{A}\f$, the preconditioner
111 * includes the interior blocks as well. The full block preconditioner
112 * is defined as:
113 *
114 * \f[\mathbf{M}^{-1}=\left[\begin{array}{cccc}
115 * \mathrm{Diag}[(\mathbf{A})_{vv}] & & & \\
116 * & (\mathbf{A})_{eb} & & \\
117 * & & (\mathbf{A})_{fb} & \\
118 * & & & (\mathbf{A})_{int} \end{array}\right] \f]
119 *
120 * where the interior portion is added at the end following the
121 * ordering inside data structures.
122 */
124{
125 ExpListSharedPtr expList = m_linsys.lock()->GetLocMat().lock();
128 DNekScalMatSharedPtr bnd_mat;
129
130 int i, j, k, n, cnt, gId;
131 int meshVertId, meshEdgeId, meshFaceId;
132
133 auto asmMap = m_locToGloMap.lock();
134
135 const int nExp = expList->GetExpSize();
136 const int nDirBnd = asmMap->GetNumGlobalDirBndCoeffs();
137
138 // Grab periodic geometry information.
139 PeriodicMap periodicVerts, periodicEdges, periodicFaces;
140 expList->GetPeriodicEntities(periodicVerts, periodicEdges, periodicFaces);
141
142 // The vectors below are of size 3(+1) to have separate storage for
143 // vertices, edges and faces.
144 // For the full matrix approach additional storage for the
145 // interior matrices is required
146 int nStorage = m_isFull ? 4 : 3;
147
148 // Maps from geometry ID to the matrix representing the extracted
149 // portion of S_1. For example idMats[2] folds the S_1 face blocks.
150 vector<map<int, vector<NekDouble>>> idMats(nStorage);
151
152 // Maps from the global ID, as obtained from AssemblyMapCG's
153 // localToGlobalMap, to the geometry ID.
154 vector<map<int, int>> gidMeshIds(3);
155
156 // Maps from the global ID to the number of degrees of freedom for
157 // this geometry object.
158 vector<map<int, int>> gidDofs(nStorage);
159
160 // Array containing maximum information needed for the universal
161 // numbering later. For i = 0,1,2 for each geometry dimension:
162 //
163 // maxVertIds[2*i] = maximum geometry ID at dimension i
164 // maxVertIds[2*i+1] = maximum number of degrees of freedom for all
165 // elements of dimension i.
166 Array<OneD, int> maxVertIds(6, -1);
167
168 // Figure out mapping from each elemental contribution to offset in
169 // (vert,edge,face) triples.
170 for (cnt = n = 0; n < nExp; ++n)
171 {
172 exp = expList->GetExp(n);
173
174 // Grab reference to
175 // StaticCond: Local Schur complement matrix
176 // Full: Boundary and interior matrix
177 DNekScalMatSharedPtr schurMat;
178
179 if (m_isFull)
180 {
181 // Get local matrix
182 auto tmpMat = m_linsys.lock()->GetBlock(n);
183
184 // Get number of boundary dofs
185 int nBndDofs = exp->NumBndryCoeffs();
186 int nIntDofs = exp->GetNcoeffs() - nBndDofs;
187
188 // Get boundary & interior maps
189 Array<OneD, unsigned int> bndMap(nBndDofs), intMap(nIntDofs);
190 exp->GetBoundaryMap(bndMap);
191 exp->GetInteriorMap(intMap);
192
193 // Create temporary matrices
194 DNekMatSharedPtr bndryMat =
196
197 // Extract boundary and send to StaticCond (Schur)
198 // framework
199 for (int i = 0; i < nBndDofs; ++i)
200 {
201 for (int j = 0; j < nBndDofs; ++j)
202 {
203 (*bndryMat)(i, j) = (*tmpMat)(bndMap[i], bndMap[j]);
204 }
205 }
206
207 schurMat =
209
210 // Extract interior matrix and save for block matrix setting
211 vector<NekDouble> tmpStore(nIntDofs * nIntDofs);
212 for (int i = 0; i < nIntDofs; ++i)
213 {
214 for (int j = 0; j < nIntDofs; ++j)
215 {
216 tmpStore[j + i * nIntDofs] =
217 (*tmpMat)(intMap[i], intMap[j]);
218 }
219 }
220
221 // Save interior mats and nDofs
222 idMats[3][n] = tmpStore;
223 gidDofs[3][n] = nIntDofs;
224 // then intMat gets added to a block diagonal matrix that
225 // handles the extra interior dofs no need to do any
226 // communication for that because all interior dofs are
227 // local to each element
228 }
229 else
230 {
231 schurMat = m_linsys.lock()->GetStaticCondBlock(n)->GetBlock(0, 0);
232 }
233
234 // Process vertices to extract relevant portion of the Schur
235 // complement matrix.
236 for (i = 0; i < exp->GetNverts(); ++i)
237 {
238 meshVertId = exp->GetGeom()->GetVid(i);
239 int locId = exp->GetVertexMap(i);
240
241 // Get the global ID of this vertex.
242 gId = asmMap->GetLocalToGlobalMap(cnt + locId) - nDirBnd;
243
244 // Ignore all Dirichlet vertices.
245 if (gId < 0)
246 {
247 continue;
248 }
249
250 gidDofs[0][gId] = 1;
251
252 // Extract vertex value from Schur complement matrix.
253 NekDouble vertVal = (*schurMat)(locId, locId);
254
255 // See if we have processed this vertex from another
256 // element.
257 auto gIt = idMats[0].find(gId);
258
259 if (gIt == idMats[0].end())
260 {
261 // If not then put our 'matrix' inside idMats.
262 idMats[0][gId] = vector<NekDouble>(1, vertVal);
263 }
264 else
265 {
266 // Otherwise combine with the value that is already
267 // there (i.e. do assembly on this degree of freedom).
268 gIt->second[0] += vertVal;
269 }
270
271 // Now check to see if the vertex is periodic. If it is,
272 // then we change meshVertId to be the minimum of all the
273 // other periodic vertices, so that we don't end up
274 // duplicating the matrix in our final block matrix.
275 auto pIt = periodicVerts.find(meshVertId);
276 if (pIt != periodicVerts.end())
277 {
278 for (j = 0; j < pIt->second.size(); ++j)
279 {
280 meshVertId = min(meshVertId, pIt->second[j].id);
281 }
282 }
283
284 // Finally record the other information we need into the
285 // other matrices.
286 gidMeshIds[0][gId] = meshVertId;
287 maxVertIds[0] = max(maxVertIds[0], meshVertId);
288 maxVertIds[1] = 1;
289 }
290
291 // Process edges. This logic is mostly the same as the previous
292 // block.
293 for (i = 0; i < exp->GetGeom()->GetNumEdges(); ++i)
294 {
295 meshEdgeId = exp->GetGeom()->GetEid(i);
296
297 Array<OneD, unsigned int> bmap, bmap2;
299 StdRegions::Orientation edgeOrient = exp->GetGeom()->GetEorient(i);
300
301 // Check if this edge is periodic. We may need to flip
302 // orientation if it is.
303 auto pIt = periodicEdges.find(meshEdgeId);
304 if (pIt != periodicEdges.end())
305 {
306 pair<int, StdRegions::Orientation> idOrient =
307 DeterminePeriodicEdgeOrientId(meshEdgeId, edgeOrient,
308 pIt->second);
309 meshEdgeId = idOrient.first;
310 edgeOrient = idOrient.second;
311 }
312
313 // Grab edge interior map, and the edge inverse boundary
314 // map, so that we can extract this edge from the Schur
315 // complement matrix.
316 if (exp->GetGeom()->GetNumFaces()) // 3D Element calls
317 {
319 ->GetEdgeInteriorToElementMap(i, bmap, sign, edgeOrient);
320 bmap2 = exp->as<LocalRegions::Expansion3D>()
321 ->GetEdgeInverseBoundaryMap(i);
322 }
323 else
324 {
325 exp->GetTraceInteriorToElementMap(i, bmap, sign, edgeOrient);
326 bmap2 = exp->as<LocalRegions::Expansion2D>()
327 ->GetTraceInverseBoundaryMap(i);
328 }
329
330 // Allocate temporary storage for the extracted edge matrix.
331 const int nEdgeCoeffs = bmap.size();
332 vector<NekDouble> tmpStore(nEdgeCoeffs * nEdgeCoeffs);
333
334 gId = asmMap->GetLocalToGlobalMap(cnt + bmap[0]);
335
336 for (j = 0; j < nEdgeCoeffs; ++j)
337 {
338 // We record the minimum ID from the edge for our
339 // maps. This follows the logic that the assembly map
340 // ordering will always give us a contiguous ordering of
341 // global degrees of freedom for edge interior
342 // coefficients.
343 gId = min(gId,
344 asmMap->GetLocalToGlobalMap(cnt + bmap[j]) - nDirBnd);
345
346 // Ignore Dirichlet edges.
347 if (gId < 0)
348 {
349 continue;
350 }
351
352 const NekDouble sign1 = sign[j];
353
354 // Extract this edge, along with sign array for assembly
355 // later.
356 for (k = 0; k < nEdgeCoeffs; ++k)
357 {
358 tmpStore[k + j * nEdgeCoeffs] =
359 sign1 * sign[k] * (*schurMat)(bmap2[j], bmap2[k]);
360 }
361 }
362
363 if (gId < 0)
364 {
365 continue;
366 }
367
368 gidDofs[1][gId] = nEdgeCoeffs;
369
370 // Assemble this edge matrix with another one, if it exists.
371 auto gIt = idMats[1].find(gId);
372 if (gIt == idMats[1].end())
373 {
374 idMats[1][gId] = tmpStore;
375 }
376 else
377 {
378 ASSERTL1(tmpStore.size() == gIt->second.size(),
379 "Number of modes mismatch");
380 Vmath::Vadd(nEdgeCoeffs * nEdgeCoeffs, &gIt->second[0], 1,
381 &tmpStore[0], 1, &gIt->second[0], 1);
382 }
383
384 gidMeshIds[1][gId] = meshEdgeId;
385 maxVertIds[2] = max(maxVertIds[2], meshEdgeId);
386 maxVertIds[3] = max(maxVertIds[3], nEdgeCoeffs);
387 }
388
389 // Process faces. This logic is mostly the same as the previous
390 // block.
391 for (i = 0; i < exp->GetGeom()->GetNumFaces(); ++i)
392 {
393 meshFaceId = exp->GetGeom()->GetFid(i);
394
395 Array<OneD, unsigned int> bmap, bmap2;
397 StdRegions::Orientation faceOrient = exp->GetGeom()->GetForient(i);
398
399 // Check if this face is periodic. We may need to flip
400 // orientation if it is.
401 auto pIt = periodicFaces.find(meshFaceId);
402 if (pIt != periodicFaces.end())
403 {
404 meshFaceId = min(meshFaceId, pIt->second[0].id);
405 faceOrient = DeterminePeriodicFaceOrient(faceOrient,
406 pIt->second[0].orient);
407 }
408
409 exp->GetTraceInteriorToElementMap(i, bmap, sign, faceOrient);
410 bmap2 = exp->as<LocalRegions::Expansion3D>()
411 ->GetTraceInverseBoundaryMap(i);
412
413 // Allocate temporary storage for the extracted face matrix.
414 const int nFaceCoeffs = bmap.size();
415 vector<NekDouble> tmpStore(nFaceCoeffs * nFaceCoeffs);
416
417 gId = asmMap->GetLocalToGlobalMap(cnt + bmap[0]);
418
419 for (j = 0; j < nFaceCoeffs; ++j)
420 {
421 gId = min(gId,
422 asmMap->GetLocalToGlobalMap(cnt + bmap[j]) - nDirBnd);
423
424 // Ignore Dirichlet faces.
425 if (gId < 0)
426 {
427 continue;
428 }
429
430 const NekDouble sign1 = sign[j];
431
432 // Extract this face, along with sign array for assembly
433 // later.
434 for (k = 0; k < nFaceCoeffs; ++k)
435 {
436 tmpStore[k + j * nFaceCoeffs] =
437 sign1 * sign[k] * (*schurMat)(bmap2[j], bmap2[k]);
438 }
439 }
440
441 if (gId < 0)
442 {
443 continue;
444 }
445
446 gidDofs[2][gId] = nFaceCoeffs;
447
448 // Assemble this face matrix with another one, if it exists.
449 auto gIt = idMats[2].find(gId);
450 if (gIt == idMats[2].end())
451 {
452 idMats[2][gId] = tmpStore;
453 }
454 else
455 {
456 ASSERTL1(tmpStore.size() == gIt->second.size(),
457 "Number of modes mismatch");
458 Vmath::Vadd(nFaceCoeffs * nFaceCoeffs, &gIt->second[0], 1,
459 &tmpStore[0], 1, &gIt->second[0], 1);
460 }
461
462 gidMeshIds[2][gId] = meshFaceId;
463 maxVertIds[4] = max(maxVertIds[4], meshFaceId);
464 maxVertIds[5] = max(maxVertIds[5], nFaceCoeffs);
465 }
466
467 cnt += exp->GetNcoeffs();
468 }
469
470 // Perform a reduction to find maximum vertex, edge and face
471 // geometry IDs.
472 m_comm = expList->GetSession()->GetComm()->GetRowComm();
473 m_comm->AllReduce(maxVertIds, LibUtilities::ReduceMax);
474
475 // Concatenate all matrices into contiguous storage and figure out
476 // universal ID numbering.
477 vector<NekDouble> storageBuf;
478 vector<long> globalToUniversal;
479
480 for (i = 0, cnt = 1; i < 3; ++i)
481 {
482 const int maxDofs = maxVertIds[2 * i + 1];
483
484 // Note that iterating over the map uses the property that keys
485 // are accessed in order of ascending order, putting everything
486 // in the right order for the global system.
487 for (auto &gIt : idMats[i])
488 {
489 // Copy matrix into storage.
490 storageBuf.insert(storageBuf.end(), gIt.second.begin(),
491 gIt.second.end());
492
493 // Get mesh ID from global ID number.
494 ASSERTL1(gidMeshIds[i].count(gIt.first) > 0,
495 "Unable to find global ID " +
496 boost::lexical_cast<string>(gIt.first) +
497 " inside map");
498 meshVertId = gidMeshIds[i][gIt.first];
499
500 for (j = 0; j < gIt.second.size(); ++j)
501 {
502 globalToUniversal.push_back(cnt +
503 meshVertId * maxDofs * maxDofs + j);
504 }
505
506 // Free up the temporary storage.
507 gIt.second.clear();
508 }
509
510 cnt += (maxVertIds[2 * i] + 1) * maxDofs * maxDofs;
511 }
512
513 ASSERTL1(storageBuf.size() == globalToUniversal.size(),
514 "Storage buffer and global to universal map size does "
515 "not match");
516
517 Array<OneD, NekDouble> storageData(storageBuf.size(), &storageBuf[0]);
518 Array<OneD, long> globalToUniversalMap(globalToUniversal.size(),
519 &globalToUniversal[0]);
520
521 // Use GS to assemble data between processors.
522 Gs::gs_data *tmpGs =
523 Gs::Init(globalToUniversalMap, m_comm,
524 expList->GetSession()->DefinesCmdLineArgument("verbose"));
525 Gs::Gather(storageData, Gs::gs_add, tmpGs);
526
527 // Figure out what storage we need in the block matrix.
528 int nblksSize = 1 + idMats[1].size() + idMats[2].size();
529 if (m_isFull)
530 {
531 nblksSize += idMats[3].size();
532 }
533 Array<OneD, unsigned int> n_blks(nblksSize);
534
535 // Vertex block is a diagonal matrix.
536 n_blks[0] = idMats[0].size();
537
538 // Now extract number of rows in each edge and face block from the
539 // gidDofs map.
540 cnt = 1;
541 for (i = 1; i < nStorage; ++i)
542 {
543 for (auto &gIt : idMats[i])
544 {
545 ASSERTL1(gidDofs[i].count(gIt.first) > 0,
546 "Unable to find number of degrees of freedom for "
547 "global ID " +
548 boost::lexical_cast<string>(gIt.first));
549
550 n_blks[cnt++] = gidDofs[i][gIt.first];
551 }
552 }
553
554 // Allocate storage for the block matrix.
555 m_blkMat =
557
558 // We deal with the vertex matrix separately since all vertices can
559 // be condensed into a single, block-diagonal matrix.
561 n_blks[0], n_blks[0], 0.0, eDIAGONAL);
562
563 // Fill the vertex matrix with the inverse of each vertex value.
564 cnt = 0;
565 for (auto gIt = idMats[0].begin(); gIt != idMats[0].end(); ++gIt, ++cnt)
566 {
567 (*vertMat)(cnt, cnt) = 1.0 / storageData[cnt];
568 }
569
570 // Put the vertex matrix in the block matrix.
571 m_blkMat->SetBlock(0, 0, vertMat);
572
573 // Finally, grab the matrices from the block storage, invert them
574 // and place them in the correct position inside the block matrix.
575 int cnt2 = 1;
576 for (i = 1; i < nStorage; ++i)
577 {
578 for (auto &gIt : idMats[i])
579 {
580 int nDofs = gidDofs[i][gIt.first];
581
582 DNekMatSharedPtr tmp =
584
585 for (j = 0; j < nDofs; ++j)
586 {
587 for (k = 0; k < nDofs; ++k)
588 {
589 // Interior matrices do not need mapping
590 // and are stored in idMats[3]
591 if (m_isFull && i == 3)
592 {
593 (*tmp)(j, k) = gIt.second[k + j * nDofs];
594 }
595 // Boundary matrices
596 else
597 {
598 (*tmp)(j, k) = storageData[k + j * nDofs + cnt];
599 }
600 }
601 }
602
603 cnt += nDofs * nDofs;
604
605 tmp->Invert();
606 m_blkMat->SetBlock(cnt2, cnt2, tmp);
607 ++cnt2;
608 }
609 }
610}
611
612/**
613 * @brief Construct a block preconditioner for the hybridized
614 * discontinuous Galerkin system.
615 *
616 * This system is built in a similar fashion to its continuous variant
617 * found in PreconditionerBlock::BlockPreconditionerCG. In this setting
618 * however, the matrix is constructed as:
619 *
620 * \f[ M^{-1} = \mathrm{Diag}[ (\mathbf{S_{1}})_{f}^{-1} ] \f]
621 *
622 * where each matrix is the Schur complement system restricted to a
623 * single face of the trace system.
624 */
626{
627 std::shared_ptr<MultiRegions::ExpList> expList =
628 ((m_linsys.lock())->GetLocMat()).lock();
629 std::shared_ptr<MultiRegions::ExpList> trace = expList->GetTrace();
632 DNekScalMatSharedPtr bnd_mat;
633
635 std::dynamic_pointer_cast<AssemblyMapDG>(m_locToGloMap.lock());
636
637 int i, j, k, n, cnt, cnt2;
638
639 // Figure out number of Dirichlet trace elements
640 int nTrace = expList->GetTrace()->GetExpSize();
641 int nDir = asmMap->GetNumGlobalDirBndCoeffs();
642
643 for (cnt = n = 0; n < nTrace; ++n)
644 {
645 if (cnt >= nDir)
646 {
647 break;
648 }
649
650 cnt += trace->GetExp(n)->GetNcoeffs();
651 }
652
653 nDir = n;
654
655 // Allocate storage for block matrix. Need number of unique faces in
656 // trace space.
657 int maxTraceSize = -1;
658 map<int, int> arrayOffsets;
659
660 for (cnt = 0, n = nDir; n < nTrace; ++n)
661 {
662 int nCoeffs = trace->GetExp(n)->GetNcoeffs();
663 int nCoeffs2 = nCoeffs * nCoeffs;
664
665 arrayOffsets[n] = cnt;
666 cnt += nCoeffs2;
667
668 if (maxTraceSize < nCoeffs2)
669 {
670 maxTraceSize = nCoeffs2;
671 }
672 }
673
674 // Find maximum trace size.
675 m_comm = expList->GetSession()->GetComm()->GetRowComm();
676 m_comm->AllReduce(maxTraceSize, LibUtilities::ReduceMax);
677
678 // Zero matrix storage.
679 Array<OneD, NekDouble> tmpStore(cnt, 0.0);
680
681 // Assemble block matrices for each trace element.
682 for (cnt = n = 0; n < expList->GetExpSize(); ++n)
683 {
684 int elmt = n;
685 locExpansion = expList->GetExp(elmt);
686
688 asmMap->GetElmtToTrace()[elmt];
689
690 // Block matrix (lambda)
691 loc_mat = (m_linsys.lock())->GetStaticCondBlock(n);
692 bnd_mat = loc_mat->GetBlock(0, 0);
693
694 int nFacets = locExpansion->GetNtraces();
695
696 for (cnt2 = i = 0; i < nFacets; ++i)
697 {
698 int nCoeffs = elmtToTraceMap[i]->GetNcoeffs();
699 int elmtId = elmtToTraceMap[i]->GetElmtId();
700
701 // Ignore Dirichlet edges.
702 if (elmtId < nDir)
703 {
704 cnt += nCoeffs;
705 cnt2 += nCoeffs;
706 continue;
707 }
708
709 NekDouble *off = &tmpStore[arrayOffsets[elmtId]];
710
711 for (j = 0; j < nCoeffs; ++j)
712 {
713 NekDouble sign1 = asmMap->GetLocalToGlobalBndSign(cnt + j);
714 for (k = 0; k < nCoeffs; ++k)
715 {
716 NekDouble sign2 = asmMap->GetLocalToGlobalBndSign(cnt + k);
717 off[j * nCoeffs + k] +=
718 (*bnd_mat)(cnt2 + j, cnt2 + k) * sign1 * sign2;
719 }
720 }
721
722 cnt += nCoeffs;
723 cnt2 += nCoeffs;
724 }
725 }
726
727 // Set up IDs for universal numbering.
728 Array<OneD, long> uniIds(tmpStore.size());
729 for (cnt = 0, n = nDir; n < nTrace; ++n)
730 {
731 LocalRegions::ExpansionSharedPtr traceExp = trace->GetExp(n);
732 int nCoeffs = traceExp->GetNcoeffs();
733 int geomId = traceExp->GetGeom()->GetGlobalID();
734
735 for (i = 0; i < nCoeffs * nCoeffs; ++i)
736 {
737 uniIds[cnt++] = geomId * maxTraceSize + i + 1;
738 }
739 }
740
741 // Assemble matrices across partitions.
742 Gs::gs_data *gsh =
743 Gs::Init(uniIds, m_comm,
744 expList->GetSession()->DefinesCmdLineArgument("verbose"));
745 Gs::Gather(tmpStore, Gs::gs_add, gsh);
746
747 // Set up diagonal block matrix
748 Array<OneD, unsigned int> n_blks(nTrace - nDir);
749 for (n = 0; n < nTrace - nDir; ++n)
750 {
751 n_blks[n] = trace->GetExp(n + nDir)->GetNcoeffs();
752 }
753
754 m_blkMat =
756
757 for (n = 0; n < nTrace - nDir; ++n)
758 {
759 int nCoeffs = trace->GetExp(n + nDir)->GetNcoeffs();
760 DNekMatSharedPtr tmp =
762 NekDouble *off = &tmpStore[arrayOffsets[n + nDir]];
763
764 for (i = 0; i < nCoeffs; ++i)
765 {
766 for (j = 0; j < nCoeffs; ++j)
767 {
768 (*tmp)(i, j) = off[i * nCoeffs + j];
769 }
770 }
771
772 tmp->Invert();
773 m_blkMat->SetBlock(n, n, tmp);
774 }
775}
776
777/**
778 * @brief Apply preconditioner to \p pInput and store the result in \p
779 * pOutput.
780 */
782 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput,
783 const bool &isLocal)
784{
785 // Get assembly map and solver type
786 auto asmMap = m_locToGloMap.lock();
787
788 // Determine matrix size
789 int nGlobal = m_isFull ? asmMap->GetNumGlobalCoeffs()
790 : asmMap->GetNumGlobalBndCoeffs();
791 int nDir = asmMap->GetNumGlobalDirBndCoeffs();
792 int nNonDir = nGlobal - nDir;
793
794 // Apply preconditioner
795 DNekBlkMat &M = (*m_blkMat);
796
797 if (isLocal)
798 {
799 Array<OneD, NekDouble> wk(nGlobal);
800 Array<OneD, NekDouble> wk1(nGlobal);
801 asmMap->Assemble(pInput, wk);
802
803 NekVector<NekDouble> r(nNonDir, wk + nDir, eWrapper);
804 NekVector<NekDouble> z(nNonDir, wk1 + nDir, eWrapper);
805
806 z = M * r;
807 Vmath::Zero(nDir, wk1, 1);
808
809 asmMap->GlobalToLocal(wk1, pOutput);
810 }
811 else
812 {
813 NekVector<NekDouble> r(nNonDir, pInput, eWrapper);
814 NekVector<NekDouble> z(nNonDir, pOutput, eWrapper);
815 z = M * r;
816 }
817}
818} // namespace Nektar::MultiRegions
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:47
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
void BlockPreconditionerCG(void)
Construct a block preconditioner from for the continuous Galerkin system.
static PreconditionerSharedPtr create(const std::shared_ptr< GlobalLinSys > &plinsys, const std::shared_ptr< AssemblyMap > &pLocToGloMap)
Creates an instance of this class.
void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &isLocal=false) override
Apply preconditioner to pInput and store the result in pOutput.
void BlockPreconditionerHDG(void)
Construct a block preconditioner for the hybridized discontinuous Galerkin system.
PreconditionerBlock(const std::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
static std::string className
Name of class.
LibUtilities::CommSharedPtr m_comm
const std::weak_ptr< GlobalLinSys > m_linsys
std::weak_ptr< AssemblyMap > m_locToGloMap
static void Gather(Nektar::Array< OneD, NekDouble > pU, gs_op pOp, gs_data *pGsh, Nektar::Array< OneD, NekDouble > pBuffer=NullNekDouble1DArray)
Performs a gather-scatter operation of the provided values.
Definition: GsLib.hpp:265
static gs_data * Init(const Nektar::Array< OneD, long > pId, const LibUtilities::CommSharedPtr &pComm, bool verbose=true)
Initialise Gather-Scatter map.
Definition: GsLib.hpp:190
@ gs_add
Definition: GsLib.hpp:60
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:66
std::shared_ptr< AssemblyMapDG > AssemblyMapDGSharedPtr
Definition: AssemblyMapDG.h:46
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
PreconFactory & GetPreconFactory()
pair< int, StdRegions::Orientation > DeterminePeriodicEdgeOrientId(int meshEdgeId, StdRegions::Orientation edgeOrient, const vector< PeriodicEntity > &periodicEdges)
Determine orientation of an edge to its periodic equivalents, as well as the ID of the representative...
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:50
std::map< int, std::vector< PeriodicEntity > > PeriodicMap
StdRegions::Orientation DeterminePeriodicFaceOrient(StdRegions::Orientation faceOrient, StdRegions::Orientation perFaceOrient)
Determine relative orientation between two faces.
std::vector< double > z(NPUPPER)
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:79
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273