21#include <opm/simulators/linalg/ParallelOverlappingILU0.hpp>
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>
31#include <opm/simulators/linalg/GraphColoring.hpp>
32#include <opm/simulators/linalg/matrixblock.hh>
43void ghost_last_bilu0_decomposition (M& A, std::size_t interiorSize)
46 assert(interiorSize <= A.N());
47 using rowiterator =
typename M::RowIterator;
48 using coliterator =
typename M::ColIterator;
49 using block =
typename M::block_type;
52 for (rowiterator i = A.begin(); i.index() < interiorSize; ++i)
55 coliterator endij=(*i).end();
59 for (ij=(*i).begin(); ij.index()<i.index(); ++ij)
62 coliterator jj = A[ij.index()].find(ij.index());
65 (*ij).rightmultiply(*jj);
68 coliterator endjk=A[ij.index()].end();
69 coliterator jk=jj; ++jk;
70 coliterator ik=ij; ++ik;
71 while (ik!=endij && jk!=endjk)
72 if (ik.index()==jk.index())
81 if (ik.index()<jk.index())
89 if (ij.index()!=i.index())
90 DUNE_THROW(Dune::ISTLError,
"diagonal entry missing");
94 catch (Dune::FMatrixError & e) {
95 DUNE_THROW(Dune::ISTLError,
"ILU failed to invert matrix block");
101template<
class M,
class CRS,
class InvVector>
102void convertToCRS(
const M& A, CRS& lower, CRS& upper, InvVector& inv)
104 OPM_TIMEBLOCK(convertToCRS);
112 using size_type =
typename M :: size_type;
117 lower.resize( A.N() );
118 upper.resize( A.N() );
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) {
131 numLower += numLowerRow;
132 numUpper += (*i).size() - numLowerRow - 1;
134 assert(numLower + numUpper + A.N() == A.nonzeroes());
136 lower.reserveAdditional( numLower );
140 size_type colcount = 0;
141 lower.rows_[ 0 ] = colcount;
142 for (
auto i=A.begin(); i!=endi; ++i, ++row)
144 const size_type iIndex = i.index();
147 for (
auto j=(*i).begin(); j.index() < iIndex; ++j )
149 lower.push_back( (*j), j.index() );
152 lower.rows_[ iIndex+1 ] = colcount;
155 assert(colcount == numLower);
157 const auto rendi = A.beforeBegin();
160 upper.rows_[ 0 ] = colcount ;
162 upper.reserveAdditional( numUpper );
166 for (
auto i=A.beforeEnd(); i!=rendi; --i, ++ row )
168 const size_type iIndex = i.index();
172 for (
auto j=(*i).beforeEnd(); j.index()>=iIndex; --j )
174 const size_type jIndex = j.index();
175 if( j.index() == iIndex )
180 else if ( j.index() >= i.index() )
182 upper.push_back( (*j), jIndex );
186 upper.rows_[ row+1 ] = colcount;
188 assert(colcount == numUpper);
192size_t set_interiorSize( [[maybe_unused]]
size_t N,
size_t interiorSize, [[maybe_unused]]
const PI& comm)
199size_t set_interiorSize(
size_t N,
size_t interiorSize,
const Dune::OwnerOverlapCopyCommunication<int,int>& comm)
203 auto indexSet = comm.indexSet();
206 for (
auto idx = indexSet.begin(); idx!=indexSet.end(); ++idx)
208 if (idx->local().attribute()==1)
210 auto loc = idx->local().local();
223template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
224Dune::SolverCategory::Category
225ParallelOverlappingILU0<Matrix,Domain,Range,ParallelInfoT>::category()
const
227 return std::is_same_v<ParallelInfoT, Dune::Amg::SequentialInformation> ?
228 Dune::SolverCategory::sequential : Dune::SolverCategory::overlapping;
231template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
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)
245 interiorSize_ = A.N();
251template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
254 const ParallelInfo& comm,
const int n,
const field_type 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)
265 interiorSize_ = A.N();
271template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
279template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
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)
293 interiorSize_ = A.N();
299template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
302 const ParallelInfo& comm,
304 size_type interiorSize,
bool redblack,
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)
317 assert(interiorSize <= A_->N());
321template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
323apply (Domain& v,
const Range& d)
325 OPM_TIMEBLOCK(
apply);
330 using dblock =
typename Range ::block_type;
331 using vblock =
typename Domain::block_type;
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())
339 OPM_THROW(std::logic_error,
"ILU: number of lower and upper rows must be the same");
343 for (size_type i = 0; i < lowerLoopEnd; ++i)
345 dblock rhs( md[ i ] );
346 const size_type rowI = lower_.rows_[ i ];
347 const size_type rowINext = lower_.rows_[ i+1 ];
349 for (size_type col = rowI; col < rowINext; ++col)
351 lower_.values_[ col ].mmv( mv[ lower_.cols_[ col ] ], rhs );
357 for (size_type i = upperLoopStart; i < iEnd; ++i)
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 ];
364 for (size_type col = rowI; col < rowINext; ++col)
366 upper_.values_[ col ].mmv( mv[ upper_.cols_[ col ] ], rhs );
370 inv_[ i ].mv( rhs, vBlock);
373 copyOwnerToAll( mv );
381template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
387 comm_->copyOwnerToAll(v, v);
391template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
395 OPM_TIMEBLOCK(update);
399 if (comm_ && comm_->communicator().size() <= 0)
403 OPM_THROW(std::logic_error,
"Expected a matrix with zero rows for an invalid communicator.");
412 int ilu_setup_successful = 1;
414 const int rank = comm_ ? comm_->communicator().rank() : 0;
418 using Graph = Dune::Amg::MatrixGraph<const Matrix>;
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_ )
436 std::vector<std::size_t> inverseOrdering(ordering_.size());
437 std::size_t index = 0;
438 for (
const auto newIndex : ordering_)
440 inverseOrdering[newIndex] = index++;
445 OPM_TIMEBLOCK(iluDecomposition);
446 if (iluIteration_ == 0) {
449 interiorSize_ = detail::set_interiorSize(A_->N(), interiorSize_, *comm_);
450 assert(interiorSize_ <= A_->N());
454 if (ordering_.empty())
457 OPM_TIMEBLOCK(iluDecompositionMakeMatrix);
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) {
472 ILU_ = std::make_unique<Matrix>(*A_);
477 ILU_ = std::make_unique<Matrix>(A_->N(), A_->M(),
478 A_->nonzeroes(), Matrix::row_wise);
481 for (
auto iter = newA.createbegin(), iend = newA.createend(); iter != iend; ++iter)
483 const auto& row = (*A_)[inverseOrdering[iter.index()]];
484 for (
auto col = row.begin(), cend = row.end(); col != cend; ++col)
486 iter.insert(ordering_[col.index()]);
490 for (
auto iter = A_->begin(), iend = A_->end(); iter != iend; ++iter)
492 auto& newRow = newA[ordering_[iter.index()]];
493 for (
auto col = iter->begin(), cend = iter->end(); col != cend; ++col)
495 newRow[ordering_[col.index()]] = *col;
503 detail::milu0_decomposition ( *ILU_);
506 detail::milu0_decomposition ( *ILU_, detail::identityFunctor<typename Matrix::field_type>,
507 detail::signFunctor<typename Matrix::field_type> );
510 detail::milu0_decomposition ( *ILU_, detail::absFunctor<typename Matrix::field_type>,
511 detail::signFunctor<typename Matrix::field_type> );
514 detail::milu0_decomposition ( *ILU_, detail::identityFunctor<typename Matrix::field_type>,
515 detail::isPositiveFunctor<typename Matrix::field_type> );
518 if (interiorSize_ == A_->N())
519 Dune::ILU::blockILU0Decomposition( *ILU_ );
521 detail::ghost_last_bilu0_decomposition(*ILU_, interiorSize_);
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())
540 milun_decomposition( *A_, iluIteration_, milu_, *ILU_, *reorderer, *inverseReorderer );
543 catch (
const Dune::MatrixBlockError& error)
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;
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)
557 throw Dune::MatrixBlockError();
561 detail::convertToCRS(*ILU_, lower_, upper_, inv_);
564template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
574 return const_cast<Range&
>(d);
588template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
608template<
class Matrix,
class Domain,
class Range,
class ParallelInfoT>
612 if (!ordering_.empty())
615 for (
const auto index : ordering_)
617 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:233
Domain & reorderV(Domain &v)
Reorder V if needed and return a reference to it.
Definition ParallelOverlappingILU0_impl.hpp:590
std::vector< std::size_t > ordering_
the reordering of the unknowns
Definition ParallelOverlappingILU0.hpp:340
Domain reorderedV_
The reordered left hand side.
Definition ParallelOverlappingILU0.hpp:344
Range reorderedD_
The reordered right hand side.
Definition ParallelOverlappingILU0.hpp:342
const field_type w_
The relaxation factor to use.
Definition ParallelOverlappingILU0.hpp:348
Range & reorderD(const Range &d)
Reorder D if needed and return a reference to it.
Definition ParallelOverlappingILU0_impl.hpp:566
void apply(Domain &v, const Range &d) override
Apply the preconditoner.
Definition ParallelOverlappingILU0_impl.hpp:323
typename Domain::field_type field_type
The field type of the preconditioner.
Definition ParallelOverlappingILU0.hpp:142
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
MILU_VARIANT
Definition MILU.hpp:34
@ MILU_1
sum(dropped entries)
Definition MILU.hpp:38
@ MILU_2
sum(dropped entries)
Definition MILU.hpp:40
@ MILU_3
sum(|dropped entries|)
Definition MILU.hpp:42
@ MILU_4
sum(dropped entries)
Definition MILU.hpp:44
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