ParallelOverlappingILU0_impl.hpp
Go to the documentation of this file.
1/*
2 Copyright 2015, 2022 Dr. Blatt - HPC-Simulation-Software & Services
3 Copyright 2015 Statoil AS
4
5 This file is part of the Open Porous Media project (OPM).
6
7 OPM is free software: you can redistribute it and/or modify
8 it under the terms of the GNU General Public License as published by
9 the Free Software Foundation, either version 3 of the License, or
10 (at your option) any later version.
11
12 OPM is distributed in the hope that it will be useful,
13 but WITHOUT ANY WARRANTY; without even the implied warranty of
14 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 GNU General Public License for more details.
16
17 You should have received a copy of the GNU General Public License
18 along with OPM. If not, see <http://www.gnu.org/licenses/>.
19*/
20
22
23#include <dune/common/version.hh>
24
25#include <dune/istl/ilu.hh>
26#include <dune/istl/owneroverlapcopy.hh>
27
28#include <opm/common/ErrorMacros.hpp>
29#include <opm/common/TimingMacros.hpp>
30
33
34#include <cassert>
35
36namespace Opm
37{
38namespace detail
39{
40
42template<class M>
43void ghost_last_bilu0_decomposition (M& A, std::size_t interiorSize)
44{
45 OPM_TIMEBLOCK(GhostLastBlockILU0Decomp);
46 // iterator types
47 assert(interiorSize <= A.N());
48 using rowiterator = typename M::RowIterator;
49 using coliterator = typename M::ColIterator;
50 using block = typename M::block_type;
51
52 // implement left looking variant with stored inverse
53 for (rowiterator i = A.begin(); i.index() < interiorSize; ++i)
54 {
55 // coliterator is diagonal after the following loop
56 coliterator endij=(*i).end(); // end of row i
57 coliterator ij;
58
59 // eliminate entries left of diagonal; store L factor
60 for (ij=(*i).begin(); ij.index()<i.index(); ++ij)
61 {
62 // find A_jj which eliminates A_ij
63 coliterator jj = A[ij.index()].find(ij.index());
64
65 // compute L_ij = A_jj^-1 * A_ij
66 (*ij).rightmultiply(*jj);
67
68 // modify row
69 coliterator endjk=A[ij.index()].end(); // end of row j
70 coliterator jk=jj; ++jk;
71 coliterator ik=ij; ++ik;
72 while (ik!=endij && jk!=endjk)
73 if (ik.index()==jk.index())
74 {
75 block B(*jk);
76 B.leftmultiply(*ij);
77 *ik -= B;
78 ++ik; ++jk;
79 }
80 else
81 {
82 if (ik.index()<jk.index())
83 ++ik;
84 else
85 ++jk;
86 }
87 }
88
89 // invert pivot and store it in A
90 if (ij.index()!=i.index())
91 DUNE_THROW(Dune::ISTLError,"diagonal entry missing");
92 try {
93 (*ij).invert(); // compute inverse of diagonal block
94 }
95 catch (Dune::FMatrixError & e) {
96 DUNE_THROW(Dune::ISTLError,"ILU failed to invert matrix block");
97 }
98 }
99}
100
102template<class M, class CRS, class InvVector>
103void convertToCRS(const M& A, CRS& lower, CRS& upper, InvVector& inv)
104{
105 OPM_TIMEBLOCK(convertToCRS);
106 // No need to do anything for 0 rows. Return to prevent indexing a
107 // a zero sized array.
108 if ( A.N() == 0 )
109 {
110 return;
111 }
112
113 using size_type = typename M :: size_type;
114
115 lower.clear();
116 upper.clear();
117 inv.clear();
118 lower.resize( A.N() );
119 upper.resize( A.N() );
120 inv.resize( A.N() );
121
122 // Count the lower and upper matrix entries.
123 size_type numLower = 0;
124 size_type numUpper = 0;
125 const auto endi = A.end();
126 for (auto i = A.begin(); i != endi; ++i) {
127 const size_type iIndex = i.index();
128 size_type numLowerRow = 0;
129 for (auto j = (*i).begin(); j.index() < iIndex; ++j) {
130 ++numLowerRow;
131 }
132 numLower += numLowerRow;
133 numUpper += (*i).size() - numLowerRow - 1;
134 }
135 assert(numLower + numUpper + A.N() == A.nonzeroes());
136
137 lower.reserveAdditional( numLower );
138
139 // implement left looking variant with stored inverse
140 size_type row = 0;
141 size_type colcount = 0;
142 lower.rows_[ 0 ] = colcount;
143 for (auto i=A.begin(); i!=endi; ++i, ++row)
144 {
145 const size_type iIndex = i.index();
146
147 // eliminate entries left of diagonal; store L factor
148 for (auto j=(*i).begin(); j.index() < iIndex; ++j )
149 {
150 lower.push_back( (*j), j.index() );
151 ++colcount;
152 }
153 lower.rows_[ iIndex+1 ] = colcount;
154 }
155
156 assert(colcount == numLower);
157
158 const auto rbegini = std::make_reverse_iterator(A.begin());
159 row = 0;
160 colcount = 0;
161 upper.rows_[ 0 ] = colcount ;
162
163 upper.reserveAdditional( numUpper );
164
165 // NOTE: upper and inv store entries in reverse order, reverse here
166 // relative to ILU
167 auto rindex = [](auto it) { return std::prev(it.base()).index(); };
168 for (auto i=std::make_reverse_iterator(A.end()); i!=rbegini; ++i, ++ row )
169 {
170 const size_type iIndex = rindex(i);
171
172 // store in reverse row order
173 // eliminate entries left of diagonal; store L factor
174 for (auto j=std::make_reverse_iterator(i->end()); rindex(j)>=iIndex; ++j )
175 {
176 const size_type jIndex = rindex(j);
177 if( rindex(j) == iIndex )
178 {
179 inv[ row ] = (*j);
180 break;
181 }
182 else if ( rindex(j) >= rindex(i) )
183 {
184 upper.push_back( (*j), jIndex );
185 ++colcount ;
186 }
187 }
188 upper.rows_[ row+1 ] = colcount;
189 }
190 assert(colcount == numUpper);
191}
192
193template <class PI>
194size_t set_interiorSize( [[maybe_unused]] size_t N, size_t interiorSize, [[maybe_unused]] const PI& comm)
195{
196 return interiorSize;
197}
198
199#if HAVE_MPI
200template<>
201size_t set_interiorSize(size_t N, size_t interiorSize, const Dune::OwnerOverlapCopyCommunication<int,int>& comm)
202{
203 if (interiorSize<N)
204 return interiorSize;
205 auto indexSet = comm.indexSet();
206
207 size_t new_is = 0;
208 for (auto idx = indexSet.begin(); idx!=indexSet.end(); ++idx)
209 {
210 if (idx->local().attribute()==1)
211 {
212 auto loc = idx->local().local();
213 if (loc > new_is) {
214 new_is = loc;
215 }
216 }
217 }
218 return new_is + 1;
219}
220#endif
221
222} // end namespace detail
223
224
225template<class Matrix, class Domain, class Range, class ParallelInfoT>
226Dune::SolverCategory::Category
228{
229 return std::is_same_v<ParallelInfoT, Dune::Amg::SequentialInformation> ?
230 Dune::SolverCategory::sequential : Dune::SolverCategory::overlapping;
231}
232
233template<class Matrix, class Domain, class Range, class ParallelInfoT>
235ParallelOverlappingILU0(const Matrix& A,
236 const int n, const field_type w,
237 MILU_VARIANT milu, bool redblack,
238 bool reorder_sphere)
239 : lower_(),
240 upper_(),
241 inv_(),
242 comm_(nullptr), w_(w),
243 relaxation_( std::abs( w - 1.0 ) > 1e-15 ),
244 A_(&reinterpret_cast<const Matrix&>(A)), iluIteration_(n),
245 milu_(milu), redBlack_(redblack), reorderSphere_(reorder_sphere)
246{
247 interiorSize_ = A.N();
248 // BlockMatrix is a Subclass of FieldMatrix that just adds
249 // methods. Therefore this cast should be safe.
250 update();
251}
252
253template<class Matrix, class Domain, class Range, class ParallelInfoT>
255ParallelOverlappingILU0(const Matrix& A,
256 const ParallelInfo& comm, const int n, const field_type w,
257 MILU_VARIANT milu, bool redblack,
258 bool reorder_sphere)
259 : lower_(),
260 upper_(),
261 inv_(),
262 comm_(&comm), w_(w),
263 relaxation_( std::abs( w - 1.0 ) > 1e-15 ),
264 A_(&reinterpret_cast<const Matrix&>(A)), iluIteration_(n),
265 milu_(milu), redBlack_(redblack), reorderSphere_(reorder_sphere)
266{
267 interiorSize_ = A.N();
268 // BlockMatrix is a Subclass of FieldMatrix that just adds
269 // methods. Therefore this cast should be safe.
270 update();
271}
272
273template<class Matrix, class Domain, class Range, class ParallelInfoT>
275ParallelOverlappingILU0(const Matrix& A,
276 const field_type w, MILU_VARIANT milu, bool redblack,
277 bool reorder_sphere)
278 : ParallelOverlappingILU0( A, 0, w, milu, redblack, reorder_sphere )
279{}
280
281template<class Matrix, class Domain, class Range, class ParallelInfoT>
283ParallelOverlappingILU0(const Matrix& A,
284 const ParallelInfo& comm, const field_type w,
285 MILU_VARIANT milu, bool redblack,
286 bool reorder_sphere)
287 : lower_(),
288 upper_(),
289 inv_(),
290 comm_(&comm), w_(w),
291 relaxation_( std::abs( w - 1.0 ) > 1e-15 ),
292 A_(&reinterpret_cast<const Matrix&>(A)), iluIteration_(0),
293 milu_(milu), redBlack_(redblack), reorderSphere_(reorder_sphere)
294{
295 interiorSize_ = A.N();
296 // BlockMatrix is a Subclass of FieldMatrix that just adds
297 // methods. Therefore this cast should be safe.
298 update();
299}
300
301template<class Matrix, class Domain, class Range, class ParallelInfoT>
303ParallelOverlappingILU0(const Matrix& A,
304 const ParallelInfo& comm,
305 const field_type w, MILU_VARIANT milu,
306 size_type interiorSize, bool redblack,
307 bool reorder_sphere)
308 : lower_(),
309 upper_(),
310 inv_(),
311 comm_(&comm), w_(w),
312 relaxation_( std::abs( w - 1.0 ) > 1e-15 ),
313 interiorSize_(interiorSize),
314 A_(&reinterpret_cast<const Matrix&>(A)), iluIteration_(0),
315 milu_(milu), redBlack_(redblack), reorderSphere_(reorder_sphere)
316{
317 // BlockMatrix is a Subclass of FieldMatrix that just adds
318 // methods. Therefore this cast should be safe.
319 assert(interiorSize <= A_->N());
320 update( );
321}
322
323template<class Matrix, class Domain, class Range, class ParallelInfoT>
325apply (Domain& v, const Range& d)
326{
327 OPM_TIMEBLOCK(apply);
328 Range& md = reorderD(d);
329 Domain& mv = reorderV(v);
330
331 // iterator types
332 using dblock = typename Range ::block_type;
333 using vblock = typename Domain::block_type;
334
335 const size_type iEnd = lower_.rows();
336 const size_type lastRow = iEnd - 1;
337 size_type upperLoopStart = iEnd - interiorSize_;
338 size_type lowerLoopEnd = interiorSize_;
339 if (iEnd != upper_.rows())
340 {
341 OPM_THROW(std::logic_error,"ILU: number of lower and upper rows must be the same");
342 }
343
344 // lower triangular solve
345 for (size_type i = 0; i < lowerLoopEnd; ++i)
346 {
347 dblock rhs( md[ i ] );
348 const size_type rowI = lower_.rows_[ i ];
349 const size_type rowINext = lower_.rows_[ i+1 ];
350
351 for (size_type col = rowI; col < rowINext; ++col)
352 {
353 lower_.values_[ col ].mmv( mv[ lower_.cols_[ col ] ], rhs );
354 }
355
356 mv[ i ] = rhs; // Lii = I
357 }
358
359 for (size_type i = upperLoopStart; i < iEnd; ++i)
360 {
361 vblock& vBlock = mv[ lastRow - i ];
362 vblock rhs ( vBlock );
363 const size_type rowI = upper_.rows_[ i ];
364 const size_type rowINext = upper_.rows_[ i+1 ];
365
366 for (size_type col = rowI; col < rowINext; ++col)
367 {
368 upper_.values_[ col ].mmv( mv[ upper_.cols_[ col ] ], rhs );
369 }
370
371 // apply inverse and store result
372 inv_[ i ].mv( rhs, vBlock);
373 }
374
375 copyOwnerToAll( mv );
376
377 if( relaxation_ ) {
378 mv *= w_;
379 }
380 reorderBack(mv, v);
381}
382
383template<class Matrix, class Domain, class Range, class ParallelInfoT>
384template<class V>
386copyOwnerToAll(V& v) const
387{
388 if( comm_ ) {
389 comm_->copyOwnerToAll(v, v);
390 }
391}
392
393template<class Matrix, class Domain, class Range, class ParallelInfoT>
395update()
396{
397 OPM_TIMEBLOCK(update);
398 // (For older DUNE versions the communicator might be
399 // invalid if redistribution in AMG happened on the coarset level.
400 // Therefore we check for nonzero size
401 if (comm_ && comm_->communicator().size() <= 0)
402 {
403 if (A_->N() > 0)
404 {
405 OPM_THROW(std::logic_error, "Expected a matrix with zero rows for an invalid communicator.");
406 }
407 else
408 {
409 // simply set the communicator to null
410 comm_ = nullptr;
411 }
412 }
413
414 int ilu_setup_successful = 1;
415 std::string message;
416 const int rank = comm_ ? comm_->communicator().rank() : 0;
417
418 if (redBlack_)
419 {
420 using Graph = Dune::Amg::MatrixGraph<const Matrix>;
421 Graph graph(*A_);
422 auto colorsTuple = colorVerticesWelshPowell(graph);
423 const auto& colors = std::get<0>(colorsTuple);
424 const auto& verticesPerColor = std::get<2>(colorsTuple);
425 auto noColors = std::get<1>(colorsTuple);
426 if ( reorderSphere_ )
427 {
428 ordering_ = reorderVerticesSpheres(colors, noColors, verticesPerColor,
429 graph, 0);
430 }
431 else
432 {
433 ordering_ = reorderVerticesPreserving(colors, noColors, verticesPerColor,
434 graph);
435 }
436 }
437
438 std::vector<std::size_t> inverseOrdering(ordering_.size());
439 {
440 OPM_TIMEBLOCK(createInverseOrdering);
441 std::size_t index = 0;
442 for (const auto newIndex : ordering_)
443 {
444 inverseOrdering[newIndex] = index++;
445 }
446 }
447
448 try
449 {
450 OPM_TIMEBLOCK(iluDecomposition);
451 if (iluIteration_ == 0) {
452
453 if (comm_) {
454 interiorSize_ = detail::set_interiorSize(A_->N(), interiorSize_, *comm_);
455 assert(interiorSize_ <= A_->N());
456 }
457
458 // create ILU-0 decomposition
459 if (ordering_.empty())
460 {
461 OPM_TIMEBLOCK(iluDecompositionUpdateMatrix);
462 if (ILU_) {
463 OPM_TIMEBLOCK(iluDecompositionCopyEntries);
464 // The ILU_ matrix is already a copy with the same
465 // sparse structure as A_, but the values of A_ may
466 // have changed, so we must copy all elements.
467 for (std::size_t row = 0; row < A_->N(); ++row) {
468 const auto& Arow = (*A_)[row];
469 auto Ait = Arow.begin();
470 auto Iit = (*ILU_)[row].begin();
471 for (; Ait != Arow.end(); ++Ait, ++Iit) {
472 *Iit = *Ait;
473 }
474 }
475 } else {
476 OPM_TIMEBLOCK(iluDecompositionDuplicateMatrix);
477 // First call, must duplicate matrix.
478 ILU_ = std::make_unique<Matrix>(*A_);
479 }
480 }
481 else
482 {
483 ILU_ = std::make_unique<Matrix>(A_->N(), A_->M(),
484 A_->nonzeroes(), Matrix::row_wise);
485 auto& newA = *ILU_;
486 // Create sparsity pattern
487 auto endcreateA = newA.createend();
488 for (auto iter = newA.createbegin(); iter != endcreateA; ++iter)
489 {
490 const auto& row = (*A_)[inverseOrdering[iter.index()]];
491 for (auto col = row.begin(), cend = row.end(); col != cend; ++col)
492 {
493 iter.insert(ordering_[col.index()]);
494 }
495 }
496 // Copy values
497 for (auto iter = A_->begin(); iter != A_->end(); ++iter)
498 {
499 auto newRow = newA.begin() + ordering_[iter.index()];
500 for (auto&& [A_ij, j] : sparseRange(*newRow))
501 {
502 (*newRow)[ordering_[j]] = A_ij;
503 }
504 }
505 }
506
507 switch (milu_)
508 {
511 break;
513 detail::milu0_decomposition ( *ILU_, detail::identityFunctor<typename Matrix::field_type>,
514 detail::signFunctor<typename Matrix::field_type> );
515 break;
517 detail::milu0_decomposition ( *ILU_, detail::absFunctor<typename Matrix::field_type>,
518 detail::signFunctor<typename Matrix::field_type> );
519 break;
521 detail::milu0_decomposition ( *ILU_, detail::identityFunctor<typename Matrix::field_type>,
522 detail::isPositiveFunctor<typename Matrix::field_type> );
523 break;
524 default:
525 if (interiorSize_ == A_->N())
526 Dune::ILU::blockILU0Decomposition( *ILU_ );
527 else
528 detail::ghost_last_bilu0_decomposition(*ILU_, interiorSize_);
529 break;
530 }
531 }
532 else {
533 // create ILU-n decomposition
534 ILU_ = std::make_unique<Matrix>(A_->N(), A_->M(), Matrix::row_wise);
535 std::unique_ptr<detail::Reorderer> reorderer, inverseReorderer;
536 if (ordering_.empty())
537 {
538 reorderer.reset(new detail::NoReorderer());
539 inverseReorderer.reset(new detail::NoReorderer());
540 }
541 else
542 {
543 reorderer.reset(new detail::RealReorderer(ordering_));
544 inverseReorderer.reset(new detail::RealReorderer(inverseOrdering));
545 }
546
547 milun_decomposition( *A_, iluIteration_, milu_, *ILU_, *reorderer, *inverseReorderer );
548 }
549 }
550 catch (const Dune::MatrixBlockError& error)
551 {
552 message = error.what();
553 std::cerr << "Exception occurred on process " << rank << " during " <<
554 "setup of ILU0 preconditioner with message: "
555 << message<<std::endl;
556 ilu_setup_successful = 0;
557 }
558
559 // Check whether there was a problem on some process
560 {
561 OPM_TIMEBLOCK(checkParallelFailure);
562 const bool parallel_failure = comm_ && comm_->communicator().min(ilu_setup_successful) == 0;
563 const bool local_failure = ilu_setup_successful == 0;
564 if (local_failure || parallel_failure)
565 {
566 throw Dune::MatrixBlockError();
567 }
568 }
569
570 // store ILU in simple CRS format
571 detail::convertToCRS(*ILU_, lower_, upper_, inv_);
572}
573
574template<class Matrix, class Domain, class Range, class ParallelInfoT>
576reorderD(const Range& d)
577{
578 if (ordering_.empty())
579 {
580 // As d is non-const in the apply method of the
581 // solver casting away constness in this particular
582 // setting is not undefined. It is ugly though but due
583 // to the preconditioner interface of dune-istl.
584 return const_cast<Range&>(d);
585 }
586 else
587 {
588 reorderedD_.resize(d.size());
589 std::size_t i = 0;
590 for (const auto index : ordering_)
591 {
592 reorderedD_[index] = d[i++];
593 }
594 return reorderedD_;
595 }
596}
597
598template<class Matrix, class Domain, class Range, class ParallelInfoT>
600reorderV(Domain& v)
601{
602 if (ordering_.empty())
603 {
604 return v;
605 }
606 else
607 {
608 reorderedV_.resize(v.size());
609 std::size_t i = 0;
610 for (const auto index : ordering_)
611 {
612 reorderedV_[index] = v[i++];
613 }
614 return reorderedV_;
615 }
616}
617
618template<class Matrix, class Domain, class Range, class ParallelInfoT>
620reorderBack(const Range& reorderedV, Range& v)
621{
622 if (!ordering_.empty())
623 {
624 std::size_t i = 0;
625 for (const auto index : ordering_)
626 {
627 v[i++] = reorderedV[index];
628 }
629 }
630}
631
632} // end namespace Opm
A two-step version of an overlapping Schwarz preconditioner using one step ILU0 as.
Definition: ParallelOverlappingILU0.hpp:131
ParallelOverlappingILU0(const Matrix &A, const int n, const field_type w, MILU_VARIANT milu, bool redblack=false, bool reorder_sphere=true)
Constructor.
Definition: ParallelOverlappingILU0_impl.hpp:235
Domain & reorderV(Domain &v)
Reorder V if needed and return a reference to it.
Definition: ParallelOverlappingILU0_impl.hpp:600
size_type interiorSize_
Definition: ParallelOverlappingILU0.hpp:350
void reorderBack(const Range &reorderedV, Range &v)
Definition: ParallelOverlappingILU0_impl.hpp:620
Range & reorderD(const Range &d)
Reorder D if needed and return a reference to it.
Definition: ParallelOverlappingILU0_impl.hpp:576
void copyOwnerToAll(V &v) const
Definition: ParallelOverlappingILU0_impl.hpp:386
Dune::SolverCategory::Category category() const override
Definition: ParallelOverlappingILU0_impl.hpp:227
void update() override
Definition: ParallelOverlappingILU0_impl.hpp:395
void apply(Domain &v, const Range &d) override
Apply the preconditoner.
Definition: ParallelOverlappingILU0_impl.hpp:325
typename matrix_type::size_type size_type
Definition: ParallelOverlappingILU0.hpp:145
typename Domain::field_type field_type
The field type of the preconditioner.
Definition: ParallelOverlappingILU0.hpp:142
void milu0_decomposition(M &A, FieldFunct< M > absFunctor=signFunctor< typename M::field_type >, FieldFunct< M > signFunctor=oneFunctor< typename M::field_type >, std::vector< typename M::block_type > *diagonal=nullptr)
void convertToCRS(const M &A, CRS &lower, CRS &upper, InvVector &inv)
compute ILU decomposition of A. A is overwritten by its decomposition
Definition: ParallelOverlappingILU0_impl.hpp:103
void milun_decomposition(const M &A, int n, MILU_VARIANT milu, M &ILU, Reorderer &ordering, Reorderer &inverseOrdering)
size_t set_interiorSize(size_t N, size_t interiorSize, const PI &comm)
Definition: ParallelOverlappingILU0_impl.hpp:194
void ghost_last_bilu0_decomposition(M &A, std::size_t interiorSize)
Compute Blocked ILU0 decomposition, when we know junk ghost rows are located at the end of A.
Definition: ParallelOverlappingILU0_impl.hpp:43
Definition: blackoilbioeffectsmodules.hh:45
MILU_VARIANT
Definition: MILU.hpp:34
@ MILU_1
sum(dropped entries)
@ MILU_2
sum(dropped entries)
@ MILU_3
sum(|dropped entries|)
@ MILU_4
sum(dropped entries)
std::vector< std::size_t > reorderVerticesPreserving(const std::vector< int > &colors, int noColors, const std::vector< std::size_t > &verticesPerColor, const Graph &graph)
! Reorder colored graph preserving order of vertices with the same color.
Definition: GraphColoring.hpp:169
std::vector< std::size_t > reorderVerticesSpheres(const std::vector< int > &colors, int noColors, const std::vector< std::size_t > &verticesPerColor, const Graph &graph, typename Graph::VertexDescriptor root)
! Reorder Vetrices in spheres
Definition: GraphColoring.hpp:189
std::tuple< std::vector< int >, int, std::vector< std::size_t > > colorVerticesWelshPowell(const Graph &graph)
Color the vertices of graph.
Definition: GraphColoring.hpp:113
Definition: MILU.hpp:76
Definition: MILU.hpp:84