3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
parallelhelpers.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_LINEAR_PARALLELHELPERS_HH
26#define DUMUX_LINEAR_PARALLELHELPERS_HH
27
28#if HAVE_MPI
29
30#include <dune/common/version.hh>
31#include <dune/geometry/dimension.hh>
32#include <dune/grid/common/datahandleif.hh>
33#include <dune/grid/common/partitionset.hh>
34#include <dune/istl/owneroverlapcopy.hh>
35#include <dune/istl/paamg/pinfo.hh>
36
37namespace Dumux {
38
44// operator that resets result to zero at constrained DOFS
45template<class LinearSolverTraits>
47{
48 using GridView = typename LinearSolverTraits::GridView;
49 using DofMapper = typename LinearSolverTraits::DofMapper;
50 enum { dofCodim = LinearSolverTraits::dofCodim };
51
52 // TODO: this is some large number (replace by limits?)
53 static constexpr std::size_t ghostMarker_ = 1<<24;
54
55 class BaseGatherScatter
56 {
57 public:
58 BaseGatherScatter(const DofMapper& mapper)
59 : mapper_(mapper) {}
60
61 template<class EntityType>
62 int index(const EntityType& e) const
63 { return mapper_.index(e); }
64
65 bool contains(int dim, int codim) const
66 { return dofCodim == codim; }
67
68#if DUNE_VERSION_GTE(DUNE_GRID,2,7)
70 bool fixedSize(int dim, int codim) const
71 { return true; }
72#else
74 bool fixedsize(int dim, int codim) const
75 { return true; }
76#endif
77
78 template<class EntityType>
79 std::size_t size(EntityType& e) const
80 { return 1; }
81
82 template<class EntityType>
83 bool isNeitherInteriorNorBorderEntity(EntityType& e) const
84 { return e.partitionType() != Dune::InteriorEntity && e.partitionType() != Dune::BorderEntity; }
85
86 private:
87 const DofMapper& mapper_;
88 };
89
93 template<class V>
94 class ConsistencyBoxGatherScatter
95 : public BaseGatherScatter,
96 public Dune::CommDataHandleIF<ConsistencyBoxGatherScatter<V>,typename V::block_type>
97 {
98 public:
99 using DataType = typename V::block_type;
100 using BaseGatherScatter::contains;
101#if DUNE_VERSION_GTE(DUNE_GRID,2,7)
102 using BaseGatherScatter::fixedSize;
103#else
104 using BaseGatherScatter::fixedsize;
105#endif
106 using BaseGatherScatter::size;
107
108 ConsistencyBoxGatherScatter(V& container, const DofMapper& mapper)
109 : BaseGatherScatter(mapper), container_(container)
110 {}
111
112 template<class MessageBuffer, class EntityType>
113 void gather(MessageBuffer& buff, const EntityType& e) const
114 {
115 buff.write(container_[this->index(e)]);
116 }
117
118 template<class MessageBuffer, class EntityType>
119 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
120 {
121 typename V::block_type block;
122 buff.read(block);
123 container_[this->index(e)] += block;
124 }
125 private:
126 V& container_;
127 };
128
129
136 class GhostGatherScatter
137 : public BaseGatherScatter,
138 public Dune::CommDataHandleIF<GhostGatherScatter,std::size_t>
139 {
140 public:
141 using DataType = std::size_t;
142 using BaseGatherScatter::contains;
143#if DUNE_VERSION_GTE(DUNE_GRID,2,7)
144 using BaseGatherScatter::fixedSize;
145#else
146 using BaseGatherScatter::fixedsize;
147#endif
148 using BaseGatherScatter::size;
149
150 GhostGatherScatter(std::vector<std::size_t>& ranks, const DofMapper& mapper)
151 : BaseGatherScatter(mapper), ranks_(ranks)
152 {}
153
154 template<class MessageBuffer, class EntityType>
155 void gather(MessageBuffer& buff, const EntityType& e) const
156 {
157 auto& data = ranks_[this->index(e)];
158 if (this->isNeitherInteriorNorBorderEntity(e))
159 data = ghostMarker_;
160 buff.write(data);
161 }
162
163 template<class MessageBuffer, class EntityType>
164 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
165 {
166 std::size_t x;
167 buff.read(x);
168 auto& data = ranks_[this->index(e)];
169 if (this->isNeitherInteriorNorBorderEntity(e))
170 data = ghostMarker_;
171 }
172 private:
173 std::vector<std::size_t>& ranks_;
174 };
175
182 class InteriorBorderGatherScatter
183 : public BaseGatherScatter,
184 public Dune::CommDataHandleIF<InteriorBorderGatherScatter,std::size_t>
185 {
186 public:
187 using DataType = std::size_t;
188 using BaseGatherScatter::contains;
189#if DUNE_VERSION_GTE(DUNE_GRID,2,7)
190 using BaseGatherScatter::fixedSize;
191#else
192 using BaseGatherScatter::fixedsize;
193#endif
194 using BaseGatherScatter::size;
195
196 InteriorBorderGatherScatter(std::vector<std::size_t>& ranks, const DofMapper& mapper)
197 : BaseGatherScatter(mapper), ranks_(ranks)
198 {}
199
200 template<class MessageBuffer, class EntityType>
201 void gather(MessageBuffer& buff, const EntityType& e) const
202 {
203 auto& data = ranks_[this->index(e)];
204 if (this->isNeitherInteriorNorBorderEntity(e))
205 data = ghostMarker_;
206 buff.write(data);
207 }
208
209 template<class MessageBuffer, class EntityType>
210 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
211 {
212 std::size_t x;
213 buff.read(x);
214 auto& data = ranks_[this->index(e)];
215 using std::min;
216 data = this->isNeitherInteriorNorBorderEntity(e) ? x : min(data,x);
217 }
218 private:
219 std::vector<std::size_t>& ranks_;
220 };
221
226 struct NeighbourGatherScatter
227 : public BaseGatherScatter,
228 public Dune::CommDataHandleIF<NeighbourGatherScatter,int>
229 {
230 using DataType = int;
231 using BaseGatherScatter::contains;
232#if DUNE_VERSION_GTE(DUNE_GRID,2,7)
233 using BaseGatherScatter::fixedSize;
234#else
235 using BaseGatherScatter::fixedsize;
236#endif
237 using BaseGatherScatter::size;
238
239 NeighbourGatherScatter(const DofMapper& mapper, int rank, std::set<int>& neighbours)
240 : BaseGatherScatter(mapper), rank_(rank), neighbours_(neighbours)
241 {}
242
243 template<class MessageBuffer, class EntityType>
244 void gather(MessageBuffer& buff, const EntityType& e) const
245 {
246 buff.write(rank_);
247 }
248
249 template<class MessageBuffer, class EntityType>
250 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
251 {
252 int x;
253 buff.read(x);
254 neighbours_.insert(x);
255 }
256 private:
257 int rank_;
258 std::set<int>& neighbours_;
259 };
260
261
266 struct SharedGatherScatter
267 : public BaseGatherScatter,
268 public Dune::CommDataHandleIF<SharedGatherScatter,int>
269 {
270 using DataType = int;
271 using BaseGatherScatter::contains;
272#if DUNE_VERSION_GTE(DUNE_GRID,2,7)
273 using BaseGatherScatter::fixedSize;
274#else
275 using BaseGatherScatter::fixedsize;
276#endif
277 using BaseGatherScatter::size;
278
279 SharedGatherScatter(std::vector<int>& shared, const DofMapper& mapper)
280 : BaseGatherScatter(mapper), shared_(shared)
281 {}
282
283 template<class MessageBuffer, class EntityType>
284 void gather(MessageBuffer& buff, EntityType& e) const
285 {
286 int data = true;
287 buff.write(data);
288 }
289
290 template<class MessageBuffer, class EntityType>
291 void scatter(MessageBuffer& buff, const EntityType &e, std::size_t n)
292 {
293 int x;
294 buff.read(x);
295 auto& data = shared_[this->index(e)];
296 data = data || x;
297 }
298 private:
299 std::vector<int>& shared_;
300 };
301
306 template<typename GlobalIndex>
307 struct GlobalIndexGatherScatter
308 : public BaseGatherScatter,
309 public Dune::CommDataHandleIF<GlobalIndexGatherScatter<GlobalIndex>, GlobalIndex>
310 {
311 using DataType = GlobalIndex;
312 using BaseGatherScatter::contains;
313#if DUNE_VERSION_GTE(DUNE_GRID,2,7)
314 using BaseGatherScatter::fixedSize;
315#else
316 using BaseGatherScatter::fixedsize;
317#endif
318 using BaseGatherScatter::size;
319
320 GlobalIndexGatherScatter(std::vector<GlobalIndex>& globalIndices, const DofMapper& mapper)
321 : BaseGatherScatter(mapper), globalIndices_(globalIndices)
322 {}
323
324 template<class MessageBuffer, class EntityType>
325 void gather(MessageBuffer& buff, const EntityType& e) const
326 {
327 buff.write(globalIndices_[this->index(e)]);
328 }
329
330 template<class MessageBuffer, class EntityType>
331 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
332 {
333 DataType x;
334 buff.read(x);
335 using std::min;
336 auto& data = globalIndices_[this->index(e)];
337 data = min(data, x);
338 }
339 private:
340 std::vector<GlobalIndex>& globalIndices_;
341 };
342
343public:
344
345 ParallelISTLHelper(const GridView& gridView, const DofMapper& mapper)
346 : gridView_(gridView), mapper_(mapper), initialized_(false)
347 {}
348
349 // \brief Initializes the markers for ghosts and owners with the correct size and values.
350 //
352 {
353 const auto rank = gridView_.comm().rank();
354 isOwned_.resize(mapper_.size(), rank);
355 // find out about ghosts
356 GhostGatherScatter ggs(isOwned_, mapper_);
357
358 if (gridView_.comm().size() > 1)
359 gridView_.communicate(ggs, Dune::InteriorBorder_All_Interface, Dune::ForwardCommunication);
360
361 isGhost_ = isOwned_;
362
363 // partition interior/border
364 InteriorBorderGatherScatter dh(isOwned_, mapper_);
365
366 if (gridView_.comm().size() > 1)
367 gridView_.communicate(dh, Dune::InteriorBorder_InteriorBorder_Interface, Dune::ForwardCommunication);
368
369 // convert vector into mask vector
370 for (auto& v : isOwned_)
371 v = (v == rank) ? 1 : 0;
372
373 initialized_ = true;
374 }
375
376 bool isGhost(std::size_t i) const
377 { return isGhost_[i] == ghostMarker_; }
378
379 // \brief Make a vector of the box model consistent.
380 template<class B, class A>
381 void makeNonOverlappingConsistent(Dune::BlockVector<B,A>& v)
382 {
383 ConsistencyBoxGatherScatter<Dune::BlockVector<B,A> > gs(v, mapper_);
384 if (gridView_.comm().size() > 1)
385 gridView_.communicate(gs, Dune::InteriorBorder_InteriorBorder_Interface,
386 Dune::ForwardCommunication);
387 }
388
395 template<class Comm>
396 void createParallelIndexSet(Comm& comm)
397 {
398 if (!initialized_)
399 {
400 // This is the first time this function is called.
401 // Therefore we need to initialize the marker vectors for ghosts and
402 // owned dofs
404 }
405
406 if (gridView_.comm().size() <= 1)
407 {
408 comm.remoteIndices().template rebuild<false>();
409 return;
410 }
411
412 // First find out which dofs we share with other processors
413 std::vector<int> isShared(mapper_.size(), false);
414
415 SharedGatherScatter sgs(isShared, mapper_);
416 gridView_.communicate(sgs, Dune::All_All_Interface, Dune::ForwardCommunication);
417
418 // Count shared dofs that we own
419 using GlobalIndex = typename Comm::ParallelIndexSet::GlobalIndex;
420 GlobalIndex count = 0;
421 for (std::size_t i = 0; i < isShared.size(); ++i)
422 if (isShared[i] && isOwned_[i] == 1)
423 ++count;
424
425 std::vector<GlobalIndex> counts(gridView_.comm().size());
426 gridView_.comm().allgather(&count, 1, counts.data());
427
428 // Compute start index start_p = \sum_{i=0}^{i<p} counts_i
429 const int rank = gridView_.comm().rank();
430 auto start = std::accumulate(counts.begin(), counts.begin() + rank, GlobalIndex(0));
431
432 std::vector<GlobalIndex> globalIndices(mapper_.size(), std::numeric_limits<GlobalIndex>::max());
433
434 for (std::size_t i = 0; i < globalIndices.size(); ++i)
435 {
436 if (isOwned_[i] == 1 && isShared[i])
437 {
438 globalIndices[i] = start; // GlobalIndex does not implement postfix ++
439 ++start;
440 }
441 }
442
443 // publish global indices for the shared DOFS to other processors.
444 GlobalIndexGatherScatter<GlobalIndex> gigs(globalIndices, mapper_);
445 gridView_.communicate(gigs, Dune::All_All_Interface, Dune::ForwardCommunication);
446
447 resizeIndexSet_(comm, globalIndices);
448
449 // Compute neighbours using communication
450 std::set<int> neighbours;
451 NeighbourGatherScatter ngs(mapper_, gridView_.comm().rank(), neighbours);
452 gridView_.communicate(ngs, Dune::All_All_Interface, Dune::ForwardCommunication);
453 comm.remoteIndices().setNeighbours(neighbours);
454
455 comm.remoteIndices().template rebuild<false>();
456 }
457
459 const DofMapper& dofMapper() const
460 { return mapper_; }
461
463 const GridView& gridView() const
464 { return gridView_; }
465
466private:
467 template<class Comm, class GlobalIndices>
468 void resizeIndexSet_(Comm& comm, const GlobalIndices& globalIndices) const
469 {
470 comm.indexSet().beginResize();
471
472 for (std::size_t localIndex = 0; localIndex < globalIndices.size(); ++localIndex)
473 {
474 const auto globalIndex = globalIndices[localIndex];
475 if (globalIndex != std::numeric_limits<typename GlobalIndices::value_type>::max())
476 {
477 const bool isOwned = isOwned_[localIndex] > 0;
478 const auto attr = getAttribute_(comm, isOwned, isGhost(localIndex));
479 using LocalIndex = typename Comm::ParallelIndexSet::LocalIndex;
480 comm.indexSet().add(globalIndex, LocalIndex{localIndex, attr});
481 }
482 }
483
484 comm.indexSet().endResize();
485 }
486
487 template<class Comm>
488 Dune::OwnerOverlapCopyAttributeSet::AttributeSet
489 getAttribute_(const Comm& comm, const bool isOwned, const bool isGhost) const
490 {
491 if (isOwned)
492 return Dune::OwnerOverlapCopyAttributeSet::owner;
493#if DUNE_VERSION_GTE(DUNE_ISTL, 2, 7)
494 else if (isGhost && (comm.category() == static_cast<int>(Dune::SolverCategory::nonoverlapping)) )
495#else
496 else if (isGhost && (comm.getSolverCategory() == static_cast<int>(Dune::SolverCategory::nonoverlapping)) )
497#endif
498 return Dune::OwnerOverlapCopyAttributeSet::overlap;
499 else
500 return Dune::OwnerOverlapCopyAttributeSet::copy;
501 }
502
503 const GridView gridView_;
504 const DofMapper& mapper_;
505 std::vector<std::size_t> isOwned_;
506 std::vector<std::size_t> isGhost_;
507 bool initialized_;
508
509}; // class ParallelISTLHelper
510
515template<class Matrix, class GridView, class DofMapper, int dofCodim>
517{
518 static constexpr int dim = GridView::dimension;
519 using Grid = typename GridView::Traits::Grid;
520 using BlockType = typename Matrix::block_type;
521 using IDS = typename Grid::Traits::GlobalIdSet;
522 using IdType = typename IDS::IdType;
523
545 template<class IsGhostFunc>
546 struct MatrixPatternExchange
547 : public Dune::CommDataHandleIF<MatrixPatternExchange<IsGhostFunc>, IdType>
548 {
550 using DataType = IdType;
551
552 MatrixPatternExchange(const DofMapper& mapper,
553 const std::map<IdType,int>& globalToLocal,
554 const std::map<int,IdType>& localToGlobal,
555 Matrix& A,
556 std::vector<std::set<int>>& sparsity,
557 const IsGhostFunc& isGhost)
558 : mapper_(mapper), idToIndex_(globalToLocal), indexToID_(localToGlobal)
559 , sparsity_(sparsity), A_(A), isGhost_(isGhost)
560 {
561 sparsity_.resize(A.N());
562 }
563
567 bool contains (int dim, int codim) const
568 { return (codim == dofCodim); }
569
570#if DUNE_VERSION_GTE(DUNE_GRID,2,7)
572 bool fixedSize(int dim, int codim) const
573 { return false; }
574#else
576 bool fixedsize(int dim, int codim) const
577 { return false; }
578#endif
579
583 template<class EntityType>
584 std::size_t size(EntityType& e) const
585 {
586 const auto rowIdx = mapper_.index(e);
587 std::size_t n = 0;
588 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
589 if (indexToID_.count(colIt.index()))
590 n++;
591
592 return n;
593 }
594
598 template<class MessageBuffer, class EntityType>
599 void gather(MessageBuffer& buff, const EntityType& e) const
600 {
601 const auto rowIdx = mapper_.index(e);
602 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
603 {
604 auto it = indexToID_.find(colIt.index());
605 if (it != indexToID_.end())
606 buff.write(it->second);
607 }
608 }
609
613 template<class MessageBuffer, class EntityType>
614 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
615 {
616 const auto rowIdx = mapper_.index(e);
617 for (std::size_t k = 0; k < n; k++)
618 {
619 IdType id;
620 buff.read(id);
621 // only add entries corresponding to border entities
622 const auto it = idToIndex_.find(id);
623 if (it != idToIndex_.end())
624 {
625 const auto colIdx = it->second;
626 if (!sparsity_[rowIdx].count(colIdx) && !isGhost_(colIdx))
627 sparsity_[rowIdx].insert(colIdx);
628 }
629 }
630 }
631
632 private:
633 const DofMapper& mapper_;
634 const std::map<IdType,int>& idToIndex_;
635 const std::map<int,IdType>& indexToID_;
636 std::vector<std::set<int>>& sparsity_;
637 Matrix& A_;
638 const IsGhostFunc& isGhost_;
639
640 }; // class MatrixPatternExchange
641
643 struct MatrixEntry
644 {
645 IdType first;
646 BlockType second;
647 MatrixEntry (const IdType& f, const BlockType& s) : first(f), second(s) {}
648 MatrixEntry () {}
649 };
650
652 struct MatrixEntryExchange
653 : public Dune::CommDataHandleIF<MatrixEntryExchange, MatrixEntry>
654 {
656 using DataType = MatrixEntry;
657
658 MatrixEntryExchange(const DofMapper& mapper,
659 const std::map<IdType,int>& globalToLocal,
660 const std::map<int,IdType>& localToGlobal,
661 Matrix& A)
662 : mapper_(mapper), idToIndex_(globalToLocal), indexToID_(localToGlobal), A_(A)
663 {}
664
668 bool contains(int dim, int codim) const
669 { return (codim == dofCodim); }
670
671#if DUNE_VERSION_GTE(DUNE_GRID,2,7)
673 bool fixedSize(int dim, int codim) const
674 { return false; }
675#else
677 bool fixedsize(int dim, int codim) const
678 { return false; }
679#endif
680
684 template<class EntityType>
685 std::size_t size(EntityType& e) const
686 {
687 const auto rowIdx = mapper_.index(e);
688 std::size_t n = 0;
689 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
690 if (indexToID_.count(colIt.index()))
691 n++;
692
693 return n;
694 }
695
699 template<class MessageBuffer, class EntityType>
700 void gather(MessageBuffer& buff, const EntityType& e) const
701 {
702 const auto rowIdx = mapper_.index(e);
703 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
704 {
705 auto it = indexToID_.find(colIt.index());
706 if (it != indexToID_.end())
707 buff.write(MatrixEntry(it->second,*colIt));
708 }
709 }
710
714 template<class MessageBuffer, class EntityType>
715 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
716 {
717 const auto rowIdx = mapper_.index(e);
718 for (std::size_t k = 0; k < n; k++)
719 {
720 MatrixEntry m;
721 buff.read(m);
722 // only add entries corresponding to border entities
723 auto it = idToIndex_.find(m.first);
724 if (it != idToIndex_.end())
725 if (A_[rowIdx].find(it->second) != A_[rowIdx].end())
726 A_[rowIdx][it->second] += m.second;
727 }
728 }
729
730 private:
731 const DofMapper& mapper_;
732 const std::map<IdType,int>& idToIndex_;
733 const std::map<int,IdType>& indexToID_;
734 Matrix& A_;
735
736 }; // class MatrixEntryExchange
737
738public:
739
740 ParallelMatrixHelper(const GridView& gridView, const DofMapper& mapper)
741 : gridView_(gridView), mapper_(mapper)
742 {
743 idToIndex_.clear();
744 indexToID_.clear();
745
746 for (const auto& entity : entities(gridView_, Dune::Codim<dofCodim>()))
747 {
748 if (entity.partitionType() == Dune::BorderEntity)
749 {
750 const int localIdx = mapper_.index(entity);
751 IdType dofIdxGlobal = gridView_.grid().globalIdSet().id(entity);
752
753 idToIndex_.emplace(dofIdxGlobal, localIdx);
754 indexToID_.emplace(localIdx, dofIdxGlobal);
755 }
756 }
757 }
758
764 template<class IsGhostFunc>
765 void extendMatrix(Matrix& A, const IsGhostFunc& isGhost)
766 {
767 if (gridView_.comm().size() <= 1)
768 return;
769
770 Matrix tmp(A);
771 std::size_t nnz = 0;
772 // get entries from other processes
773 std::vector<std::set<int>> sparsity;
774 MatrixPatternExchange<IsGhostFunc> datahandle(mapper_, idToIndex_, indexToID_, A, sparsity, isGhost);
775 gridView_.communicate(datahandle, Dune::InteriorBorder_InteriorBorder_Interface,
776 Dune::ForwardCommunication);
777 // add own entries, count number of nonzeros
778 for (auto rowIt = A.begin(); rowIt != A.end(); ++rowIt)
779 {
780 for (auto colIt = A[rowIt.index()].begin(); colIt != A[rowIt.index()].end(); ++colIt)
781 if (!sparsity[rowIt.index()].count(colIt.index()))
782 sparsity[rowIt.index()].insert(colIt.index());
783
784 nnz += sparsity[rowIt.index()].size();
785 }
786
787 A.setSize(tmp.N(), tmp.N(), nnz);
788 A.setBuildMode(Matrix::row_wise);
789 auto citer = A.createbegin();
790 for (auto i = sparsity.begin(), end = sparsity.end(); i!=end; ++i, ++citer)
791 {
792 for (auto si = i->begin(), send = i->end(); si!=send; ++si)
793 citer.insert(*si);
794 }
795
796 // set matrix old values
797 A = 0;
798 for (auto rowIt = tmp.begin(); rowIt != tmp.end(); ++rowIt)
799 for (auto colIt = tmp[rowIt.index()].begin(); colIt != tmp[rowIt.index()].end(); ++colIt)
800 A[rowIt.index()][colIt.index()] = tmp[rowIt.index()][colIt.index()];
801 }
802
807 void sumEntries(Matrix& A)
808 {
809 if (gridView_.comm().size() <= 1)
810 return;
811
812 MatrixEntryExchange datahandle(mapper_, idToIndex_, indexToID_, A);
813 gridView_.communicate(datahandle, Dune::InteriorBorder_InteriorBorder_Interface,
814 Dune::ForwardCommunication);
815 }
816
817private:
818 const GridView gridView_;
819 const DofMapper& mapper_;
820 std::map<IdType, int> idToIndex_;
821 std::map<int, IdType> indexToID_;
822
823}; // class EntityExchanger
824
828template<class LinearSolverTraits, class ParallelTraits,
829 class Matrix, class Vector, class ParallelHelper>
830void prepareLinearAlgebraParallel(Matrix& A, Vector& b,
831 std::shared_ptr<typename ParallelTraits::Comm>& comm,
832 std::shared_ptr<typename ParallelTraits::LinearOperator>& fop,
833 std::shared_ptr<typename ParallelTraits::ScalarProduct>& sp,
834 ParallelHelper& pHelper)
835{
836 if constexpr (ParallelTraits::isNonOverlapping)
837 {
838 // extend the matrix pattern such that it is usable for a parallel solver
839 // and make right-hand side consistent
840 using GridView = typename LinearSolverTraits::GridView;
841 using DofMapper = typename LinearSolverTraits::DofMapper;
842 static constexpr int dofCodim = LinearSolverTraits::dofCodim;
843 ParallelMatrixHelper<Matrix, GridView, DofMapper, dofCodim> matrixHelper(pHelper.gridView(), pHelper.dofMapper());
844 matrixHelper.extendMatrix(A, [&pHelper](auto idx){ return pHelper.isGhost(idx); });
845 matrixHelper.sumEntries(A);
846
847 pHelper.makeNonOverlappingConsistent(b);
848
849 // create commicator, operator, scalar product
850 const auto category = Dune::SolverCategory::nonoverlapping;
851 comm = std::make_shared<typename ParallelTraits::Comm>(pHelper.gridView().comm(), category);
852 pHelper.createParallelIndexSet(*comm);
853 fop = std::make_shared<typename ParallelTraits::LinearOperator>(A, *comm);
854 sp = std::make_shared<typename ParallelTraits::ScalarProduct>(*comm);
855 }
856 else
857 {
858 // create commicator, operator, scalar product
859 const auto category = Dune::SolverCategory::overlapping;
860 comm = std::make_shared<typename ParallelTraits::Comm>(pHelper.gridView().comm(), category);
861 pHelper.createParallelIndexSet(*comm);
862 fop = std::make_shared<typename ParallelTraits::LinearOperator>(A, *comm);
863 sp = std::make_shared<typename ParallelTraits::ScalarProduct>(*comm);
864 }
865}
866
867} // end namespace Dumux
868
869#endif // HAVE_MPI
870#endif // DUMUX_PARALLELHELPERS_HH
int start(int argc, char **argv, void(*usage)(const char *, const std::string &))
Provides a main function with error handling.
Definition: start.hh:127
Definition: adapt.hh:29
LinearSolverTraitsImpl< GridGeometry, GridGeometry::discMethod > LinearSolverTraits
The type traits required for using the IstlFactoryBackend.
Definition: linearsolvertraits.hh:72
void prepareLinearAlgebraParallel(Matrix &A, Vector &b, std::shared_ptr< typename ParallelTraits::Comm > &comm, std::shared_ptr< typename ParallelTraits::LinearOperator > &fop, std::shared_ptr< typename ParallelTraits::ScalarProduct > &sp, ParallelHelper &pHelper)
Prepare linear algebra variables for parallel solvers.
Definition: parallelhelpers.hh:830
A parallel helper class providing a nonoverlapping decomposition of all degrees of freedom.
Definition: parallelhelpers.hh:47
void makeNonOverlappingConsistent(Dune::BlockVector< B, A > &v)
Definition: parallelhelpers.hh:381
bool isGhost(std::size_t i) const
Definition: parallelhelpers.hh:376
const DofMapper & dofMapper() const
Return the dofMapper.
Definition: parallelhelpers.hh:459
void createParallelIndexSet(Comm &comm)
Creates a parallel index set.
Definition: parallelhelpers.hh:396
const GridView & gridView() const
Return the gridView.
Definition: parallelhelpers.hh:463
ParallelISTLHelper(const GridView &gridView, const DofMapper &mapper)
Definition: parallelhelpers.hh:345
void initGhostsAndOwners()
Definition: parallelhelpers.hh:351
Helper class for adding up matrix entries on border.
Definition: parallelhelpers.hh:517
ParallelMatrixHelper(const GridView &gridView, const DofMapper &mapper)
Definition: parallelhelpers.hh:740
void extendMatrix(Matrix &A, const IsGhostFunc &isGhost)
communicates values for the sparsity pattern of the new matrix.
Definition: parallelhelpers.hh:765
void sumEntries(Matrix &A)
Sums up the entries corresponding to border vertices.
Definition: parallelhelpers.hh:807