mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Work around unstable matrix inversion in DUNE 2.[34]
The versions are missing the specialized code for inverting a 3x3 matrix that makes the algorithms quite a bit more stable. With this patch we fall back to using our own MatrixBlock that does not suffer from this deficiency.
This commit is contained in:
parent
f0b60cfeed
commit
1bb0968283
@ -289,7 +289,14 @@ namespace Opm
|
||||
}
|
||||
}
|
||||
|
||||
#if DUNE_VERSION_NEWER(DUNE_ISTL, 2 , 5)
|
||||
typedef ParallelOverlappingILU0<Matrix,Vector,Vector> SeqPreconditioner;
|
||||
#else
|
||||
typedef ParallelOverlappingILU0<Dune::BCRSMatrix<Dune::MatrixBlock<typename Matrix::field_type,
|
||||
Matrix::block_type::rows,
|
||||
Matrix::block_type::cols> >,
|
||||
Vector, Vector> SeqPreconditioner;
|
||||
#endif
|
||||
|
||||
template <class Operator>
|
||||
std::unique_ptr<SeqPreconditioner> constructPrecond(Operator& opA, const Dune::Amg::SequentialInformation&) const
|
||||
@ -302,7 +309,14 @@ namespace Opm
|
||||
|
||||
#if HAVE_MPI
|
||||
typedef Dune::OwnerOverlapCopyCommunication<int, int> Comm;
|
||||
#if DUNE_VERSION_NEWER(DUNE_ISTL, 2 , 5)
|
||||
typedef ParallelOverlappingILU0<Matrix,Vector,Vector,Comm> ParPreconditioner;
|
||||
#else
|
||||
typedef ParallelOverlappingILU0<Dune::BCRSMatrix<Dune::MatrixBlock<typename Matrix::field_type,
|
||||
Matrix::block_type::rows,
|
||||
Matrix::block_type::cols> >,
|
||||
Vector, Vector, Comm> ParPreconditioner;
|
||||
#endif
|
||||
template <class Operator>
|
||||
std::unique_ptr<ParPreconditioner>
|
||||
constructPrecond(Operator& opA, const Comm& comm) const
|
||||
|
@ -238,14 +238,17 @@ public:
|
||||
\param n ILU fill in level (for testing). This does not work in parallel.
|
||||
\param w The relaxation factor.
|
||||
*/
|
||||
ParallelOverlappingILU0 (const Matrix& A, const int n, const field_type w )
|
||||
template<class Matrix1>
|
||||
ParallelOverlappingILU0 (const Matrix1& A, const int n, const field_type w )
|
||||
: lower_(),
|
||||
upper_(),
|
||||
inv_(),
|
||||
comm_(nullptr), w_(w),
|
||||
relaxation_( std::abs( w - 1.0 ) > 1e-15 )
|
||||
{
|
||||
init( A, n );
|
||||
// BlockMatrix is a Subclass of FieldMatrix that just adds
|
||||
// methods. Therefore this cast should be safe.
|
||||
init( reinterpret_cast<const Matrix&>(A), n );
|
||||
}
|
||||
|
||||
/*! \brief Constructor.
|
||||
@ -254,7 +257,8 @@ public:
|
||||
\param A The matrix to operate on.
|
||||
\param w The relaxation factor.
|
||||
*/
|
||||
ParallelOverlappingILU0 (const Matrix& A, const field_type w)
|
||||
template<class Matrix1>
|
||||
ParallelOverlappingILU0 (const Matrix1& A, const field_type w)
|
||||
: ParallelOverlappingILU0( A, 0, w )
|
||||
{
|
||||
}
|
||||
@ -266,14 +270,17 @@ public:
|
||||
\param comm communication object, e.g. Dune::OwnerOverlapCopyCommunication
|
||||
\param w The relaxation factor.
|
||||
*/
|
||||
ParallelOverlappingILU0 (const Matrix& A, const ParallelInfo& comm, const field_type w)
|
||||
template<class Matrix1>
|
||||
ParallelOverlappingILU0 (const Matrix1& A, const ParallelInfo& comm, const field_type w)
|
||||
: lower_(),
|
||||
upper_(),
|
||||
inv_(),
|
||||
comm_(&comm), w_(w),
|
||||
relaxation_( std::abs( w - 1.0 ) > 1e-15 )
|
||||
{
|
||||
init( A, 0 );
|
||||
// BlockMatrix is a Subclass of FieldMatrix that just adds
|
||||
// methods. Therefore this cast should be safe.
|
||||
init( reinterpret_cast<const Matrix&>(A), 0 );
|
||||
}
|
||||
|
||||
/*!
|
||||
|
Loading…
Reference in New Issue
Block a user