3.4
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#include <dune/common/exceptions.hh>
29#include <dune/geometry/dimension.hh>
30#include <dune/grid/common/datahandleif.hh>
31#include <dune/grid/common/partitionset.hh>
32#include <dune/istl/owneroverlapcopy.hh>
33#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
69 bool fixedSize(int dim, int codim) const
70 { return true; }
71
72 template<class EntityType>
73 std::size_t size(EntityType& e) const
74 { return 1; }
75
76 template<class EntityType>
77 bool isNeitherInteriorNorBorderEntity(EntityType& e) const
78 { return e.partitionType() != Dune::InteriorEntity && e.partitionType() != Dune::BorderEntity; }
79
80 private:
81 const DofMapper& mapper_;
82 };
83
90 class GhostGatherScatter
91 : public BaseGatherScatter,
92 public Dune::CommDataHandleIF<GhostGatherScatter,std::size_t>
93 {
94 public:
95 using DataType = std::size_t;
96 using BaseGatherScatter::contains;
97 using BaseGatherScatter::fixedSize;
98 using BaseGatherScatter::size;
99
100 GhostGatherScatter(std::vector<std::size_t>& ranks, const DofMapper& mapper)
101 : BaseGatherScatter(mapper), ranks_(ranks)
102 {}
103
104 template<class MessageBuffer, class EntityType>
105 void gather(MessageBuffer& buff, const EntityType& e) const
106 {
107 auto& data = ranks_[this->index(e)];
108 if (this->isNeitherInteriorNorBorderEntity(e))
109 data = ghostMarker_;
110 buff.write(data);
111 }
112
113 template<class MessageBuffer, class EntityType>
114 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
115 {
116 std::size_t x;
117 buff.read(x);
118 auto& data = ranks_[this->index(e)];
119 if (this->isNeitherInteriorNorBorderEntity(e))
120 data = ghostMarker_;
121 }
122 private:
123 std::vector<std::size_t>& ranks_;
124 };
125
132 class InteriorBorderGatherScatter
133 : public BaseGatherScatter,
134 public Dune::CommDataHandleIF<InteriorBorderGatherScatter,std::size_t>
135 {
136 public:
137 using DataType = std::size_t;
138 using BaseGatherScatter::contains;
139 using BaseGatherScatter::fixedSize;
140 using BaseGatherScatter::size;
141
142 InteriorBorderGatherScatter(std::vector<std::size_t>& ranks, const DofMapper& mapper)
143 : BaseGatherScatter(mapper), ranks_(ranks)
144 {}
145
146 template<class MessageBuffer, class EntityType>
147 void gather(MessageBuffer& buff, const EntityType& e) const
148 {
149 auto& data = ranks_[this->index(e)];
150 if (this->isNeitherInteriorNorBorderEntity(e))
151 data = ghostMarker_;
152 buff.write(data);
153 }
154
155 template<class MessageBuffer, class EntityType>
156 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
157 {
158 std::size_t x;
159 buff.read(x);
160 auto& data = ranks_[this->index(e)];
161 using std::min;
162 data = this->isNeitherInteriorNorBorderEntity(e) ? x : min(data,x);
163 }
164 private:
165 std::vector<std::size_t>& ranks_;
166 };
167
172 struct NeighbourGatherScatter
173 : public BaseGatherScatter,
174 public Dune::CommDataHandleIF<NeighbourGatherScatter,int>
175 {
176 using DataType = int;
177 using BaseGatherScatter::contains;
178 using BaseGatherScatter::fixedSize;
179 using BaseGatherScatter::size;
180
181 NeighbourGatherScatter(const DofMapper& mapper, int rank, std::set<int>& neighbours)
182 : BaseGatherScatter(mapper), rank_(rank), neighbours_(neighbours)
183 {}
184
185 template<class MessageBuffer, class EntityType>
186 void gather(MessageBuffer& buff, const EntityType& e) const
187 {
188 buff.write(rank_);
189 }
190
191 template<class MessageBuffer, class EntityType>
192 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
193 {
194 int x;
195 buff.read(x);
196 neighbours_.insert(x);
197 }
198 private:
199 int rank_;
200 std::set<int>& neighbours_;
201 };
202
203
208 struct SharedGatherScatter
209 : public BaseGatherScatter,
210 public Dune::CommDataHandleIF<SharedGatherScatter,int>
211 {
212 using DataType = int;
213 using BaseGatherScatter::contains;
214 using BaseGatherScatter::fixedSize;
215 using BaseGatherScatter::size;
216
217 SharedGatherScatter(std::vector<int>& shared, const DofMapper& mapper)
218 : BaseGatherScatter(mapper), shared_(shared)
219 {}
220
221 template<class MessageBuffer, class EntityType>
222 void gather(MessageBuffer& buff, EntityType& e) const
223 {
224 int data = true;
225 buff.write(data);
226 }
227
228 template<class MessageBuffer, class EntityType>
229 void scatter(MessageBuffer& buff, const EntityType &e, std::size_t n)
230 {
231 int x;
232 buff.read(x);
233 auto& data = shared_[this->index(e)];
234 data = data || x;
235 }
236 private:
237 std::vector<int>& shared_;
238 };
239
244 template<typename GlobalIndex>
245 struct GlobalIndexGatherScatter
246 : public BaseGatherScatter,
247 public Dune::CommDataHandleIF<GlobalIndexGatherScatter<GlobalIndex>, GlobalIndex>
248 {
249 using DataType = GlobalIndex;
250 using BaseGatherScatter::contains;
251 using BaseGatherScatter::fixedSize;
252 using BaseGatherScatter::size;
253
254 GlobalIndexGatherScatter(std::vector<GlobalIndex>& globalIndices, const DofMapper& mapper)
255 : BaseGatherScatter(mapper), globalIndices_(globalIndices)
256 {}
257
258 template<class MessageBuffer, class EntityType>
259 void gather(MessageBuffer& buff, const EntityType& e) const
260 {
261 buff.write(globalIndices_[this->index(e)]);
262 }
263
264 template<class MessageBuffer, class EntityType>
265 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
266 {
267 DataType x;
268 buff.read(x);
269 using std::min;
270 auto& data = globalIndices_[this->index(e)];
271 data = min(data, x);
272 }
273 private:
274 std::vector<GlobalIndex>& globalIndices_;
275 };
276
277public:
278
279 ParallelISTLHelper(const GridView& gridView, const DofMapper& mapper)
280 : gridView_(gridView), mapper_(mapper)
281 {
282 if constexpr (Detail::canCommunicate<typename GridView::Traits::Grid, dofCodim>)
283 initGhostsAndOwners_();
284 else
285 DUNE_THROW(Dune::InvalidStateException, "Cannot initialize parallel helper for a grid that cannot communicate codim-" << dofCodim << "-entities.");
286 }
287
288 bool isGhost(std::size_t i) const
289 { return isGhost_[i] == ghostMarker_; }
290
297 template<class Comm>
298 void createParallelIndexSet(Comm& comm) const
299 {
300 if constexpr (Detail::canCommunicate<typename GridView::Traits::Grid, dofCodim>)
301 {
302 if (gridView_.comm().size() <= 1)
303 {
304 comm.remoteIndices().template rebuild<false>();
305 return;
306 }
307
308 // First find out which dofs we share with other processors
309 std::vector<int> isShared(mapper_.size(), false);
310
311 SharedGatherScatter sgs(isShared, mapper_);
312 gridView_.communicate(sgs, Dune::All_All_Interface, Dune::ForwardCommunication);
313
314 // Count shared dofs that we own
315 using GlobalIndex = typename Comm::ParallelIndexSet::GlobalIndex;
316 GlobalIndex count = 0;
317 for (std::size_t i = 0; i < isShared.size(); ++i)
318 if (isShared[i] && isOwned_[i] == 1)
319 ++count;
320
321 std::vector<GlobalIndex> counts(gridView_.comm().size());
322 gridView_.comm().allgather(&count, 1, counts.data());
323
324 // Compute start index start_p = \sum_{i=0}^{i<p} counts_i
325 const int rank = gridView_.comm().rank();
326 auto start = std::accumulate(counts.begin(), counts.begin() + rank, GlobalIndex(0));
327
328 std::vector<GlobalIndex> globalIndices(mapper_.size(), std::numeric_limits<GlobalIndex>::max());
329
330 for (std::size_t i = 0; i < globalIndices.size(); ++i)
331 {
332 if (isOwned_[i] == 1 && isShared[i])
333 {
334 globalIndices[i] = start; // GlobalIndex does not implement postfix ++
335 ++start;
336 }
337 }
338
339 // publish global indices for the shared DOFS to other processors.
340 GlobalIndexGatherScatter<GlobalIndex> gigs(globalIndices, mapper_);
341 gridView_.communicate(gigs, Dune::All_All_Interface, Dune::ForwardCommunication);
342
343 resizeIndexSet_(comm, globalIndices);
344
345 // Compute neighbours using communication
346 std::set<int> neighbours;
347 NeighbourGatherScatter ngs(mapper_, gridView_.comm().rank(), neighbours);
348 gridView_.communicate(ngs, Dune::All_All_Interface, Dune::ForwardCommunication);
349 comm.remoteIndices().setNeighbours(neighbours);
350
351 comm.remoteIndices().template rebuild<false>();
352 }
353 else
354 DUNE_THROW(Dune::InvalidStateException, "Cannot build parallel index set for a grid that cannot communicate codim-" << dofCodim << "-entities.");
355 }
356
358 const DofMapper& dofMapper() const
359 { return mapper_; }
360
362 const GridView& gridView() const
363 { return gridView_; }
364
365private:
366 void initGhostsAndOwners_()
367 {
368 const auto rank = gridView_.comm().rank();
369 isOwned_.assign(mapper_.size(), rank);
370 // find out about ghosts
371 GhostGatherScatter ggs(isOwned_, mapper_);
372
373 if (gridView_.comm().size() > 1)
374 gridView_.communicate(ggs, Dune::InteriorBorder_All_Interface, Dune::ForwardCommunication);
375
376 isGhost_ = isOwned_;
377
378 // partition interior/border
379 InteriorBorderGatherScatter dh(isOwned_, mapper_);
380
381 if (gridView_.comm().size() > 1)
382 gridView_.communicate(dh, Dune::InteriorBorder_InteriorBorder_Interface, Dune::ForwardCommunication);
383
384 // convert vector into mask vector
385 for (auto& v : isOwned_)
386 v = (v == rank) ? 1 : 0;
387 }
388
389 template<class Comm, class GlobalIndices>
390 void resizeIndexSet_(Comm& comm, const GlobalIndices& globalIndices) const
391 {
392 comm.indexSet().beginResize();
393
394 for (std::size_t localIndex = 0; localIndex < globalIndices.size(); ++localIndex)
395 {
396 const auto globalIndex = globalIndices[localIndex];
397 if (globalIndex != std::numeric_limits<typename GlobalIndices::value_type>::max())
398 {
399 const bool isOwned = isOwned_[localIndex] > 0;
400 const auto attr = getAttribute_(comm, isOwned, isGhost(localIndex));
401 using LocalIndex = typename Comm::ParallelIndexSet::LocalIndex;
402 comm.indexSet().add(globalIndex, LocalIndex{localIndex, attr});
403 }
404 }
405
406 comm.indexSet().endResize();
407 }
408
409 template<class Comm>
410 Dune::OwnerOverlapCopyAttributeSet::AttributeSet
411 getAttribute_(const Comm& comm, const bool isOwned, const bool isGhost) const
412 {
413 if (isOwned)
414 return Dune::OwnerOverlapCopyAttributeSet::owner;
415 else if (isGhost && (comm.category() == static_cast<int>(Dune::SolverCategory::nonoverlapping)) )
416 return Dune::OwnerOverlapCopyAttributeSet::overlap;
417 else
418 return Dune::OwnerOverlapCopyAttributeSet::copy;
419 }
420
421 const GridView gridView_;
422 const DofMapper& mapper_;
423 std::vector<std::size_t> isOwned_;
424 std::vector<std::size_t> isGhost_;
425
426}; // class ParallelISTLHelper
427
428template<class GridView, class DofMapper, int dofCodim>
430{
431public:
432 ParallelVectorHelper(const GridView& gridView, const DofMapper& mapper)
433 : gridView_(gridView), mapper_(mapper)
434 {}
435
436 // \brief Make a vector of the box model consistent.
437 template<class Block, class Alloc>
438 void makeNonOverlappingConsistent(Dune::BlockVector<Block, Alloc>& v) const
439 {
440 if constexpr (Detail::canCommunicate<typename GridView::Traits::Grid, dofCodim>)
441 {
443 if (gridView_.comm().size() > 1)
444 gridView_.communicate(gs, Dune::InteriorBorder_InteriorBorder_Interface,
445 Dune::ForwardCommunication);
446 }
447 else
448 DUNE_THROW(Dune::InvalidStateException, "Cannot call makeNonOverlappingConsistent for a grid that cannot communicate codim-" << dofCodim << "-entities.");
449 }
450
451private:
452 const GridView gridView_;
453 const DofMapper& mapper_;
454};
455
460template<class Matrix, class GridView, class DofMapper, int dofCodim>
462{
463 static constexpr int dim = GridView::dimension;
464 using Grid = typename GridView::Traits::Grid;
465 using BlockType = typename Matrix::block_type;
466 using IDS = typename Grid::Traits::GlobalIdSet;
467 using IdType = typename IDS::IdType;
468
490 template<class IsGhostFunc>
491 struct MatrixPatternExchange
492 : public Dune::CommDataHandleIF<MatrixPatternExchange<IsGhostFunc>, IdType>
493 {
495 using DataType = IdType;
496
497 MatrixPatternExchange(const DofMapper& mapper,
498 const std::map<IdType,int>& globalToLocal,
499 const std::map<int,IdType>& localToGlobal,
500 Matrix& A,
501 std::vector<std::set<int>>& sparsity,
502 const IsGhostFunc& isGhost)
503 : mapper_(mapper), idToIndex_(globalToLocal), indexToID_(localToGlobal)
504 , sparsity_(sparsity), A_(A), isGhost_(isGhost)
505 {
506 sparsity_.resize(A.N());
507 }
508
512 bool contains (int dim, int codim) const
513 { return (codim == dofCodim); }
514
516 bool fixedSize(int dim, int codim) const
517 { return false; }
518
522 template<class EntityType>
523 std::size_t size(EntityType& e) const
524 {
525 const auto rowIdx = mapper_.index(e);
526 std::size_t n = 0;
527 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
528 if (indexToID_.count(colIt.index()))
529 n++;
530
531 return n;
532 }
533
537 template<class MessageBuffer, class EntityType>
538 void gather(MessageBuffer& buff, const EntityType& e) const
539 {
540 const auto rowIdx = mapper_.index(e);
541 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
542 {
543 auto it = indexToID_.find(colIt.index());
544 if (it != indexToID_.end())
545 buff.write(it->second);
546 }
547 }
548
552 template<class MessageBuffer, class EntityType>
553 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
554 {
555 const auto rowIdx = mapper_.index(e);
556 for (std::size_t k = 0; k < n; k++)
557 {
558 IdType id;
559 buff.read(id);
560 // only add entries corresponding to border entities
561 const auto it = idToIndex_.find(id);
562 if (it != idToIndex_.end())
563 {
564 const auto colIdx = it->second;
565 if (!sparsity_[rowIdx].count(colIdx) && !isGhost_(colIdx))
566 sparsity_[rowIdx].insert(colIdx);
567 }
568 }
569 }
570
571 private:
572 const DofMapper& mapper_;
573 const std::map<IdType,int>& idToIndex_;
574 const std::map<int,IdType>& indexToID_;
575 std::vector<std::set<int>>& sparsity_;
576 Matrix& A_;
577 const IsGhostFunc& isGhost_;
578
579 }; // class MatrixPatternExchange
580
582 struct MatrixEntry
583 {
584 IdType first;
585 BlockType second;
586 MatrixEntry (const IdType& f, const BlockType& s) : first(f), second(s) {}
587 MatrixEntry () {}
588 };
589
591 struct MatrixEntryExchange
592 : public Dune::CommDataHandleIF<MatrixEntryExchange, MatrixEntry>
593 {
595 using DataType = MatrixEntry;
596
597 MatrixEntryExchange(const DofMapper& mapper,
598 const std::map<IdType,int>& globalToLocal,
599 const std::map<int,IdType>& localToGlobal,
600 Matrix& A)
601 : mapper_(mapper), idToIndex_(globalToLocal), indexToID_(localToGlobal), A_(A)
602 {}
603
607 bool contains(int dim, int codim) const
608 { return (codim == dofCodim); }
609
611 bool fixedSize(int dim, int codim) const
612 { return false; }
613
617 template<class EntityType>
618 std::size_t size(EntityType& e) const
619 {
620 const auto rowIdx = mapper_.index(e);
621 std::size_t n = 0;
622 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
623 if (indexToID_.count(colIt.index()))
624 n++;
625
626 return n;
627 }
628
632 template<class MessageBuffer, class EntityType>
633 void gather(MessageBuffer& buff, const EntityType& e) const
634 {
635 const auto rowIdx = mapper_.index(e);
636 for (auto colIt = A_[rowIdx].begin(); colIt != A_[rowIdx].end(); ++colIt)
637 {
638 auto it = indexToID_.find(colIt.index());
639 if (it != indexToID_.end())
640 buff.write(MatrixEntry(it->second,*colIt));
641 }
642 }
643
647 template<class MessageBuffer, class EntityType>
648 void scatter(MessageBuffer& buff, const EntityType& e, std::size_t n)
649 {
650 const auto rowIdx = mapper_.index(e);
651 for (std::size_t k = 0; k < n; k++)
652 {
653 MatrixEntry m;
654 buff.read(m);
655 // only add entries corresponding to border entities
656 auto it = idToIndex_.find(m.first);
657 if (it != idToIndex_.end())
658 if (A_[rowIdx].find(it->second) != A_[rowIdx].end())
659 A_[rowIdx][it->second] += m.second;
660 }
661 }
662
663 private:
664 const DofMapper& mapper_;
665 const std::map<IdType,int>& idToIndex_;
666 const std::map<int,IdType>& indexToID_;
667 Matrix& A_;
668
669 }; // class MatrixEntryExchange
670
671public:
672
673 ParallelMatrixHelper(const GridView& gridView, const DofMapper& mapper)
674 : gridView_(gridView), mapper_(mapper)
675 {
676 idToIndex_.clear();
677 indexToID_.clear();
678
679 for (const auto& entity : entities(gridView_, Dune::Codim<dofCodim>()))
680 {
681 if (entity.partitionType() == Dune::BorderEntity)
682 {
683 const int localIdx = mapper_.index(entity);
684 IdType dofIdxGlobal = gridView_.grid().globalIdSet().id(entity);
685
686 idToIndex_.emplace(dofIdxGlobal, localIdx);
687 indexToID_.emplace(localIdx, dofIdxGlobal);
688 }
689 }
690 }
691
697 template<class IsGhostFunc>
698 void extendMatrix(Matrix& A, const IsGhostFunc& isGhost)
699 {
700 if constexpr (Detail::canCommunicate<Grid, dofCodim>)
701 {
702 if (gridView_.comm().size() <= 1)
703 return;
704
705 Matrix tmp(A);
706 std::size_t nnz = 0;
707 // get entries from other processes
708 std::vector<std::set<int>> sparsity;
709 MatrixPatternExchange<IsGhostFunc> datahandle(mapper_, idToIndex_, indexToID_, A, sparsity, isGhost);
710 gridView_.communicate(datahandle, Dune::InteriorBorder_InteriorBorder_Interface,
711 Dune::ForwardCommunication);
712 // add own entries, count number of nonzeros
713 for (auto rowIt = A.begin(); rowIt != A.end(); ++rowIt)
714 {
715 for (auto colIt = A[rowIt.index()].begin(); colIt != A[rowIt.index()].end(); ++colIt)
716 if (!sparsity[rowIt.index()].count(colIt.index()))
717 sparsity[rowIt.index()].insert(colIt.index());
718
719 nnz += sparsity[rowIt.index()].size();
720 }
721
722 A.setSize(tmp.N(), tmp.N(), nnz);
723 A.setBuildMode(Matrix::row_wise);
724 auto citer = A.createbegin();
725 for (auto i = sparsity.begin(), end = sparsity.end(); i!=end; ++i, ++citer)
726 {
727 for (auto si = i->begin(), send = i->end(); si!=send; ++si)
728 citer.insert(*si);
729 }
730
731 // set matrix old values
732 A = 0;
733 for (auto rowIt = tmp.begin(); rowIt != tmp.end(); ++rowIt)
734 for (auto colIt = tmp[rowIt.index()].begin(); colIt != tmp[rowIt.index()].end(); ++colIt)
735 A[rowIt.index()][colIt.index()] = tmp[rowIt.index()][colIt.index()];
736 }
737 else
738 DUNE_THROW(Dune::InvalidStateException, "Cannot call extendMatrix for a grid that cannot communicate codim-" << dofCodim << "-entities.");
739 }
740
745 void sumEntries(Matrix& A)
746 {
747 if constexpr (Detail::canCommunicate<Grid, dofCodim>)
748 {
749 if (gridView_.comm().size() <= 1)
750 return;
751
752 MatrixEntryExchange datahandle(mapper_, idToIndex_, indexToID_, A);
753 gridView_.communicate(datahandle, Dune::InteriorBorder_InteriorBorder_Interface,
754 Dune::ForwardCommunication);
755 }
756 else
757 DUNE_THROW(Dune::InvalidStateException, "Cannot call sumEntries for a grid that cannot communicate codim-" << dofCodim << "-entities.");
758 }
759
760private:
761 const GridView gridView_;
762 const DofMapper& mapper_;
763 std::map<IdType, int> idToIndex_;
764 std::map<int, IdType> indexToID_;
765
766};
767
771template<class LinearSolverTraits, class ParallelTraits,
772 class Matrix, class Vector, class ParallelHelper>
773void prepareLinearAlgebraParallel(Matrix& A, Vector& b, ParallelHelper& pHelper)
774{
775 if constexpr (ParallelTraits::isNonOverlapping)
776 {
777 // extend the matrix pattern such that it is usable for a parallel solver
778 // and make right-hand side consistent
779 using GridView = typename LinearSolverTraits::GridView;
780 using DofMapper = typename LinearSolverTraits::DofMapper;
781 static constexpr int dofCodim = LinearSolverTraits::dofCodim;
782 ParallelMatrixHelper<Matrix, GridView, DofMapper, dofCodim> matrixHelper(pHelper.gridView(), pHelper.dofMapper());
783 matrixHelper.extendMatrix(A, [&pHelper](auto idx){ return pHelper.isGhost(idx); });
784 matrixHelper.sumEntries(A);
785
786 ParallelVectorHelper<GridView, DofMapper, dofCodim> vectorHelper(pHelper.gridView(), pHelper.dofMapper());
787 vectorHelper.makeNonOverlappingConsistent(b);
788 }
789}
790
794template<class LinearSolverTraits, class ParallelTraits,
795 class Matrix, class Vector, class ParallelHelper>
796void prepareLinearAlgebraParallel(Matrix& A, Vector& b,
797 std::shared_ptr<typename ParallelTraits::Comm>& comm,
798 std::shared_ptr<typename ParallelTraits::LinearOperator>& fop,
799 std::shared_ptr<typename ParallelTraits::ScalarProduct>& sp,
800 ParallelHelper& pHelper)
801{
802 prepareLinearAlgebraParallel<LinearSolverTraits, ParallelTraits>(A, b, pHelper);
803 const auto category = ParallelTraits::isNonOverlapping ?
804 Dune::SolverCategory::nonoverlapping : Dune::SolverCategory::overlapping;
805
806 comm = std::make_shared<typename ParallelTraits::Comm>(pHelper.gridView().comm(), category);
807 pHelper.createParallelIndexSet(*comm);
808 fop = std::make_shared<typename ParallelTraits::LinearOperator>(A, *comm);
809 sp = std::make_shared<typename ParallelTraits::ScalarProduct>(*comm);
810}
811
812} // end namespace Dumux
813
814#endif
dune-grid capabilities compatibility layer
Contains a class to exchange entries of a vector.
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
void prepareLinearAlgebraParallel(Matrix &A, Vector &b, ParallelHelper &pHelper)
Prepare linear algebra variables for parallel solvers.
Definition: parallelhelpers.hh:773
LinearSolverTraitsImpl< GridGeometry, GridGeometry::discMethod > LinearSolverTraits
The type traits required for using the IstlFactoryBackend.
Definition: linearsolvertraits.hh:45
A parallel helper class providing a nonoverlapping decomposition of all degrees of freedom.
Definition: parallelhelpers.hh:47
bool isGhost(std::size_t i) const
Definition: parallelhelpers.hh:288
void createParallelIndexSet(Comm &comm) const
Creates a parallel index set.
Definition: parallelhelpers.hh:298
const DofMapper & dofMapper() const
Return the dofMapper.
Definition: parallelhelpers.hh:358
const GridView & gridView() const
Return the gridView.
Definition: parallelhelpers.hh:362
ParallelISTLHelper(const GridView &gridView, const DofMapper &mapper)
Definition: parallelhelpers.hh:279
Definition: parallelhelpers.hh:430
ParallelVectorHelper(const GridView &gridView, const DofMapper &mapper)
Definition: parallelhelpers.hh:432
void makeNonOverlappingConsistent(Dune::BlockVector< Block, Alloc > &v) const
Definition: parallelhelpers.hh:438
Helper class for adding up matrix entries on border.
Definition: parallelhelpers.hh:462
ParallelMatrixHelper(const GridView &gridView, const DofMapper &mapper)
Definition: parallelhelpers.hh:673
void extendMatrix(Matrix &A, const IsGhostFunc &isGhost)
communicates values for the sparsity pattern of the new matrix.
Definition: parallelhelpers.hh:698
void sumEntries(Matrix &A)
Sums up the entries corresponding to border vertices.
Definition: parallelhelpers.hh:745
A data handle class to exchange entries of a vector.
Definition: vectorcommdatahandle.hh:78