23#include <dune/common/version.hh>
25#include <dune/istl/ilu.hh>
26#include <dune/istl/owneroverlapcopy.hh>
28#include <opm/common/ErrorMacros.hpp>
29#include <opm/common/TimingMacros.hpp>
45 OPM_TIMEBLOCK(GhostLastBlockILU0Decomp);
47 assert(interiorSize <= A.N());
48 using rowiterator =
typename M::RowIterator;
49 using coliterator =
typename M::ColIterator;
50 using block =
typename M::block_type;
53 for (rowiterator i = A.begin(); i.index() < interiorSize; ++i)
56 coliterator endij=(*i).end();
60 for (ij=(*i).begin(); ij.index()<i.index(); ++ij)
63 coliterator jj = A[ij.index()].find(ij.index());
66 (*ij).rightmultiply(*jj);
69 coliterator endjk=A[ij.index()].end();
70 coliterator jk=jj; ++jk;
71 coliterator ik=ij; ++ik;
72 while (ik!=endij && jk!=endjk)
73 if (ik.index()==jk.index())
82 if (ik.index()<jk.index())
90 if (ij.index()!=i.index())
91 DUNE_THROW(Dune::ISTLError,
"diagonal entry missing");
95 catch (Dune::FMatrixError & e) {
96 DUNE_THROW(Dune::ISTLError,
"ILU failed to invert matrix block");
102template<
class M,
class CRS,
class InvVector>
113 using size_type =
typename M :: size_type;
118 lower.resize( A.N() );
119 upper.resize( A.N() );
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) {
132 numLower += numLowerRow;
133 numUpper += (*i).size() - numLowerRow - 1;
135 assert(numLower + numUpper + A.N() == A.nonzeroes());
137 lower.reserveAdditional( numLower );
141 size_type colcount = 0;
142 lower.rows_[ 0 ] = colcount;
143 for (
auto i=A.begin(); i!=endi; ++i, ++row)
145 const size_type iIndex = i.index();
148 for (
auto j=(*i).begin(); j.index() < iIndex; ++j )
150 lower.push_back( (*j), j.index() );
153 lower.rows_[ iIndex+1 ] = colcount;
156 assert(colcount == numLower);
158 const auto rendi = A.beforeBegin();
161 upper.rows_[ 0 ] = colcount ;
163 upper.reserveAdditional( numUpper );
167 for (
auto i=A.beforeEnd(); i!=rendi; --i, ++ row )
169 const size_type iIndex = i.index();
173 for (
auto j=(*i).beforeEnd(); j.index()>=iIndex; --j )
175 const size_type jIndex = j.index();
176 if( j.index() == iIndex )
181 else if ( j.index() >= i.index() )
183 upper.push_back( (*j), jIndex );
187 upper.rows_[ row+1 ] = colcount;
189 assert(colcount == numUpper);
193size_t set_interiorSize( [[maybe_unused]]
size_t N,
size_t interiorSize, [[maybe_unused]]
const PI& comm)
200size_t set_interiorSize(
size_t N,
size_t interiorSize,
const Dune::OwnerOverlapCopyCommunication<int,int>& comm)
204 auto indexSet = comm.indexSet();
207 for (
auto idx = indexSet.begin(); idx!=indexSet.end(); ++idx)
209 if (idx->local().attribute()==1)
211 auto loc = idx->local().local();
224template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
225Dune::SolverCategory::Category
228 return std::is_same_v<ParallelInfoT, Dune::Amg::SequentialInformation> ?
229 Dune::SolverCategory::sequential : Dune::SolverCategory::overlapping;
232template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
241 comm_(nullptr), w_(w),
242 relaxation_( std::abs( w - 1.0 ) > 1e-15 ),
243 A_(&reinterpret_cast<const Matrix&>(A)), iluIteration_(n),
244 milu_(milu), redBlack_(redblack), reorderSphere_(reorder_sphere)
252template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
255 const ParallelInfo& comm,
const int n,
const field_type w,
262 relaxation_( std::abs( w - 1.0 ) > 1e-15 ),
263 A_(&reinterpret_cast<const Matrix&>(A)), iluIteration_(n),
264 milu_(milu), redBlack_(redblack), reorderSphere_(reorder_sphere)
272template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
280template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
290 relaxation_( std::abs( w - 1.0 ) > 1e-15 ),
291 A_(&reinterpret_cast<const Matrix&>(A)), iluIteration_(0),
292 milu_(milu), redBlack_(redblack), reorderSphere_(reorder_sphere)
300template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
303 const ParallelInfo& comm,
311 relaxation_( std::abs( w - 1.0 ) > 1e-15 ),
312 interiorSize_(interiorSize),
313 A_(&reinterpret_cast<const Matrix&>(A)), iluIteration_(0),
314 milu_(milu), redBlack_(redblack), reorderSphere_(reorder_sphere)
318 assert(interiorSize <= A_->N());
322template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
324apply (Domain& v,
const Range& d)
326 OPM_TIMEBLOCK(apply);
327 Range& md = reorderD(d);
328 Domain& mv = reorderV(v);
331 using dblock =
typename Range ::block_type;
332 using vblock =
typename Domain::block_type;
336 size_type upperLoopStart = iEnd - interiorSize_;
338 if (iEnd != upper_.rows())
340 OPM_THROW(std::logic_error,
"ILU: number of lower and upper rows must be the same");
344 for (
size_type i = 0; i < lowerLoopEnd; ++i)
346 dblock rhs( md[ i ] );
347 const size_type rowI = lower_.rows_[ i ];
348 const size_type rowINext = lower_.rows_[ i+1 ];
350 for (
size_type col = rowI; col < rowINext; ++col)
352 lower_.values_[ col ].mmv( mv[ lower_.cols_[ col ] ], rhs );
358 for (
size_type i = upperLoopStart; i < iEnd; ++i)
360 vblock& vBlock = mv[ lastRow - i ];
361 vblock rhs ( vBlock );
362 const size_type rowI = upper_.rows_[ i ];
363 const size_type rowINext = upper_.rows_[ i+1 ];
365 for (
size_type col = rowI; col < rowINext; ++col)
367 upper_.values_[ col ].mmv( mv[ upper_.cols_[ col ] ], rhs );
371 inv_[ i ].mv( rhs, vBlock);
374 copyOwnerToAll( mv );
382template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
388 comm_->copyOwnerToAll(v, v);
392template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
396 OPM_TIMEBLOCK(update);
400 if (comm_ && comm_->communicator().size() <= 0)
404 OPM_THROW(std::logic_error,
"Expected a matrix with zero rows for an invalid communicator.");
413 int ilu_setup_successful = 1;
415 const int rank = comm_ ? comm_->communicator().rank() : 0;
419 using Graph = Dune::Amg::MatrixGraph<const Matrix>;
422 const auto& colors = std::get<0>(colorsTuple);
423 const auto& verticesPerColor = std::get<2>(colorsTuple);
424 auto noColors = std::get<1>(colorsTuple);
425 if ( reorderSphere_ )
437 std::vector<std::size_t> inverseOrdering(ordering_.size());
439 OPM_TIMEBLOCK(createInverseOrdering);
440 std::size_t index = 0;
441 for (
const auto newIndex : ordering_)
443 inverseOrdering[newIndex] = index++;
449 OPM_TIMEBLOCK(iluDecomposition);
450 if (iluIteration_ == 0) {
454 assert(interiorSize_ <= A_->N());
458 if (ordering_.empty())
460 OPM_TIMEBLOCK(iluDecompositionUpdateMatrix);
462 OPM_TIMEBLOCK(iluDecompositionCopyEntries);
466 for (std::size_t row = 0; row < A_->N(); ++row) {
467 const auto& Arow = (*A_)[row];
468 auto& ILUrow = (*ILU_)[row];
469 auto Ait = Arow.begin();
470 auto Iit = ILUrow.begin();
471 for (; Ait != Arow.end(); ++Ait, ++Iit) {
476 OPM_TIMEBLOCK(iluDecompositionDuplicateMatrix);
478 ILU_ = std::make_unique<Matrix>(*A_);
483 ILU_ = std::make_unique<Matrix>(A_->N(), A_->M(),
484 A_->nonzeroes(), Matrix::row_wise);
487 for (
auto iter = newA.createbegin(), iend = newA.createend(); iter != iend; ++iter)
489 const auto& row = (*A_)[inverseOrdering[iter.index()]];
490 for (
auto col = row.begin(), cend = row.end(); col != cend; ++col)
492 iter.insert(ordering_[col.index()]);
496 for (
auto iter = A_->begin(), iend = A_->end(); iter != iend; ++iter)
498 auto& newRow = newA[ordering_[iter.index()]];
499 for (
auto col = iter->begin(), cend = iter->end(); col != cend; ++col)
501 newRow[ordering_[col.index()]] = *col;
513 detail::signFunctor<typename Matrix::field_type> );
517 detail::signFunctor<typename Matrix::field_type> );
521 detail::isPositiveFunctor<typename Matrix::field_type> );
524 if (interiorSize_ == A_->N())
525 Dune::ILU::blockILU0Decomposition( *ILU_ );
533 ILU_ = std::make_unique<Matrix>(A_->N(), A_->M(), Matrix::row_wise);
534 std::unique_ptr<detail::Reorderer> reorderer, inverseReorderer;
535 if (ordering_.empty())
549 catch (
const Dune::MatrixBlockError& error)
551 message = error.what();
552 std::cerr <<
"Exception occurred on process " << rank <<
" during " <<
553 "setup of ILU0 preconditioner with message: "
554 << message<<std::endl;
555 ilu_setup_successful = 0;
560 OPM_TIMEBLOCK(checkParallelFailure);
561 const bool parallel_failure = comm_ && comm_->communicator().min(ilu_setup_successful) == 0;
562 const bool local_failure = ilu_setup_successful == 0;
563 if (local_failure || parallel_failure)
565 throw Dune::MatrixBlockError();
573template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
577 if (ordering_.empty())
583 return const_cast<Range&
>(d);
587 reorderedD_.resize(d.size());
589 for (
const auto index : ordering_)
591 reorderedD_[index] = d[i++];
597template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
601 if (ordering_.empty())
607 reorderedV_.resize(v.size());
609 for (
const auto index : ordering_)
611 reorderedV_[index] = v[i++];
617template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
621 if (!ordering_.empty())
624 for (
const auto index : ordering_)
626 v[i++] = reorderedV[index];
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:234
Domain & reorderV(Domain &v)
Reorder V if needed and return a reference to it.
Definition: ParallelOverlappingILU0_impl.hpp:599
size_type interiorSize_
Definition: ParallelOverlappingILU0.hpp:350
void reorderBack(const Range &reorderedV, Range &v)
Definition: ParallelOverlappingILU0_impl.hpp:619
Range & reorderD(const Range &d)
Reorder D if needed and return a reference to it.
Definition: ParallelOverlappingILU0_impl.hpp:575
void copyOwnerToAll(V &v) const
Definition: ParallelOverlappingILU0_impl.hpp:385
Dune::SolverCategory::Category category() const override
Definition: ParallelOverlappingILU0_impl.hpp:226
void update() override
Definition: ParallelOverlappingILU0_impl.hpp:394
void apply(Domain &v, const Range &d) override
Apply the preconditoner.
Definition: ParallelOverlappingILU0_impl.hpp:324
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:193
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