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