/* Copyright 2016 IRIS AS This file is part of the Open Porous Media project (OPM). OPM is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. OPM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with OPM. If not, see . */ #ifndef OPM_ISTLSOLVER_HEADER_INCLUDED #define OPM_ISTLSOLVER_HEADER_INCLUDED #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include namespace Dune { namespace ISTLUtility { //! invert matrix by calling FMatrixHelp::invert template static inline void invertMatrix (FieldMatrix &matrix) { FieldMatrix A ( matrix ); FMatrixHelp::invertMatrix(A, matrix ); } //! invert matrix by calling FMatrixHelp::invert template static inline void invertMatrix (FieldMatrix &matrix) { FieldMatrix A ( matrix ); FMatrixHelp::invertMatrix(A, matrix ); } //! invert matrix by calling FMatrixHelp::invert template static inline void invertMatrix (FieldMatrix &matrix) { FieldMatrix A ( matrix ); FMatrixHelp::invertMatrix(A, matrix ); } //! invert matrix by calling matrix.invert template static inline void invertMatrix (FieldMatrix &matrix) { matrix.invert(); } } // end ISTLUtility template class MatrixBlock : public Dune::FieldMatrix { public: typedef Dune::FieldMatrix BaseType; using BaseType :: operator= ; using BaseType :: rows; using BaseType :: cols; explicit MatrixBlock( const Scalar scalar = 0 ) : BaseType( scalar ) {} void invert() { ISTLUtility::invertMatrix( *this ); } const BaseType& asBase() const { return static_cast< const BaseType& > (*this); } BaseType& asBase() { return static_cast< BaseType& > (*this); } }; template void print_row (std::ostream& s, const MatrixBlock& A, typename FieldMatrix::size_type I, typename FieldMatrix::size_type J, typename FieldMatrix::size_type therow, int width, int precision) { print_row(s, A.asBase(), I, J, therow, width, precision); } template K& firstmatrixelement (MatrixBlock& A) { return firstmatrixelement( A.asBase() ); } template struct MatrixDimension< MatrixBlock< Scalar, n, m > > : public MatrixDimension< typename MatrixBlock< Scalar, n, m >::BaseType > { }; #if HAVE_UMFPACK /// \brief UMFPack specialization for MatrixBlock to make AMG happy /// /// Without this the empty default implementation would be used. template class UMFPack, A> > : public UMFPack, A> > { typedef UMFPack, A> > Base; typedef BCRSMatrix, A> Matrix; public: typedef BCRSMatrix, A> RealMatrix; UMFPack(const RealMatrix& matrix, int verbose, bool) : Base(reinterpret_cast(matrix), verbose) {} }; #endif #if HAVE_SUPERLU /// \brief SuperLU specialization for MatrixBlock to make AMG happy /// /// Without this the empty default implementation would be used. template class SuperLU, A> > : public SuperLU, A> > { typedef SuperLU, A> > Base; typedef BCRSMatrix, A> Matrix; public: typedef BCRSMatrix, A> RealMatrix; SuperLU(const RealMatrix& matrix, int verbose, bool reuse=true) : Base(reinterpret_cast(matrix), verbose, reuse) {} }; #endif } // end namespace Dune namespace Opm { /// This class solves the fully implicit black-oil system by /// solving the reduced system (after eliminating well variables) /// as a block-structured matrix (one block for all cell variables) for a fixed /// number of cell variables np . /// \tparam MatrixBlockType The type of the matrix block used. /// \tparam VectorBlockType The type of the vector block used. /// \tparam pressureIndex The index of the pressure component in the vector /// vector block. It is used to guide the AMG coarsening. /// Default is zero. template < class MatrixBlockType, class VectorBlockType, int pressureIndex=0 > class ISTLSolver : public NewtonIterationBlackoilInterface { typedef typename MatrixBlockType :: field_type Scalar; typedef Dune::BCRSMatrix Matrix; typedef Dune::BlockVector Vector; public: typedef Dune::AssembledLinearOperator< Matrix, Vector, Vector > AssembledLinearOperatorType; typedef NewtonIterationBlackoilInterface :: SolutionVector SolutionVector; /// Construct a system solver. /// \param[in] param parameters controlling the behaviour of the linear solvers /// \param[in] parallelInformation In the case of a parallel run /// with dune-istl the information about the parallelization. ISTLSolver(const NewtonIterationBlackoilInterleavedParameters& param, const boost::any& parallelInformation_arg=boost::any()) : iterations_( 0 ), parallelInformation_(parallelInformation_arg), isIORank_(isIORank(parallelInformation_arg)), parameters_( param ) { } /// Construct a system solver. /// \param[in] param ParameterGroup controlling the behaviour of the linear solvers /// \param[in] parallelInformation In the case of a parallel run /// with dune-istl the information about the parallelization. ISTLSolver(const ParameterGroup& param, const boost::any& parallelInformation_arg=boost::any()) : iterations_( 0 ), parallelInformation_(parallelInformation_arg), isIORank_(isIORank(parallelInformation_arg)), parameters_( param ) { } // dummy method that is not implemented for this class SolutionVector computeNewtonIncrement(const LinearisedBlackoilResidual&) const { OPM_THROW(std::logic_error,"This method is not implemented"); return SolutionVector(); } /// Solve the system of linear equations Ax = b, with A being the /// combined derivative matrix of the residual and b /// being the residual itself. /// \param[in] residual residual object containing A and b. /// \return the solution x /// \copydoc NewtonIterationBlackoilInterface::iterations int iterations () const { return iterations_; } /// \copydoc NewtonIterationBlackoilInterface::parallelInformation const boost::any& parallelInformation() const { return parallelInformation_; } public: /// \brief construct the CPR preconditioner and the solver. /// \tparam P The type of the parallel information. /// \param parallelInformation the information about the parallelization. template void constructPreconditionerAndSolve(LinearOperator& linearOperator, Vector& x, Vector& istlb, const POrComm& parallelInformation_arg, Dune::InverseOperatorResult& result) const { // Construct scalar product. typedef Dune::ScalarProductChooser ScalarProductChooser; typedef std::unique_ptr SPPointer; SPPointer sp(ScalarProductChooser::construct(parallelInformation_arg)); // Communicate if parallel. parallelInformation_arg.copyOwnerToAll(istlb, istlb); #if FLOW_SUPPORT_AMG // activate AMG if either flow_ebos is used or UMFPack is not available if( parameters_.linear_solver_use_amg_ ) { typedef ISTLUtility::CPRSelector< Matrix, Vector, Vector, POrComm> CPRSelectorType; typedef typename CPRSelectorType::AMG AMG; typedef typename CPRSelectorType::Operator MatrixOperator; std::unique_ptr< AMG > amg; std::unique_ptr< MatrixOperator > opA; if( ! std::is_same< LinearOperator, MatrixOperator > :: value ) { // create new operator in case linear operator and matrix operator differ opA.reset( CPRSelectorType::makeOperator( linearOperator.getmat(), parallelInformation_arg ) ); } const double relax = 1.0; // Construct preconditioner. constructAMGPrecond( linearOperator, parallelInformation_arg, amg, opA, relax ); // Solve. solve(linearOperator, x, istlb, *sp, *amg, result); } else #endif { // Construct preconditioner. auto precond = constructPrecond(linearOperator, parallelInformation_arg); // Solve. solve(linearOperator, x, istlb, *sp, *precond, result); } } #if DUNE_VERSION_NEWER(DUNE_ISTL, 2 , 5) typedef ParallelOverlappingILU0 SeqPreconditioner; #else typedef ParallelOverlappingILU0 >, Vector, Vector> SeqPreconditioner; #endif template std::unique_ptr constructPrecond(Operator& opA, const Dune::Amg::SequentialInformation&) const { const double relax = parameters_.ilu_relaxation_; const int ilu_fillin = parameters_.ilu_fillin_level_; std::unique_ptr precond(new SeqPreconditioner(opA.getmat(), ilu_fillin, relax)); return precond; } #if HAVE_MPI typedef Dune::OwnerOverlapCopyCommunication Comm; #if DUNE_VERSION_NEWER(DUNE_ISTL, 2 , 5) typedef ParallelOverlappingILU0 ParPreconditioner; #else typedef ParallelOverlappingILU0 >, Vector, Vector, Comm> ParPreconditioner; #endif template std::unique_ptr constructPrecond(Operator& opA, const Comm& comm) const { typedef std::unique_ptr Pointer; const double relax = parameters_.ilu_relaxation_; return Pointer(new ParPreconditioner(opA.getmat(), comm, relax)); } #endif template void constructAMGPrecond(LinearOperator& /* linearOperator */, const POrComm& comm, std::unique_ptr< AMG >& amg, std::unique_ptr< MatrixOperator >& opA, const double relax ) const { ISTLUtility::template createAMGPreconditionerPointer( *opA, relax, comm, amg ); } template void constructAMGPrecond(MatrixOperator& opA, const POrComm& comm, std::unique_ptr< AMG >& amg, std::unique_ptr< MatrixOperator >&, const double relax ) const { ISTLUtility::template createAMGPreconditionerPointer( opA, relax, comm, amg ); } /// \brief Solve the system using the given preconditioner and scalar product. template void solve(Operator& opA, Vector& x, Vector& istlb, ScalarProd& sp, Precond& precond, Dune::InverseOperatorResult& result) const { // TODO: Revise when linear solvers interface opm-core is done // Construct linear solver. // GMRes solver int verbosity = ( isIORank_ ) ? parameters_.linear_solver_verbosity_ : 0; if ( parameters_.newton_use_gmres_ ) { Dune::RestartedGMResSolver linsolve(opA, sp, precond, parameters_.linear_solver_reduction_, parameters_.linear_solver_restart_, parameters_.linear_solver_maxiter_, verbosity); // Solve system. linsolve.apply(x, istlb, result); } else { // BiCGstab solver Dune::BiCGSTABSolver linsolve(opA, sp, precond, parameters_.linear_solver_reduction_, parameters_.linear_solver_maxiter_, verbosity); // Solve system. linsolve.apply(x, istlb, result); } } /// Solve the linear system Ax = b, with A being the /// combined derivative matrix of the residual and b /// being the residual itself. /// \param[in] A matrix A /// \param[inout] x solution to be computed x /// \param[in] b right hand side b void solve(Matrix& A, Vector& x, Vector& b ) const { // Parallel version is deactivated until we figure out how to do it properly. #if HAVE_MPI if (parallelInformation_.type() == typeid(ParallelISTLInformation)) { typedef Dune::OwnerOverlapCopyCommunication Comm; const ParallelISTLInformation& info = boost::any_cast( parallelInformation_); Comm istlComm(info.communicator()); // Construct operator, scalar product and vectors needed. typedef Dune::OverlappingSchwarzOperator Operator; Operator opA(A, istlComm); solve( opA, x, b, istlComm ); } else #endif { // Construct operator, scalar product and vectors needed. Dune::MatrixAdapter< Matrix, Vector, Vector> opA( A ); solve( opA, x, b ); } } /// Solve the linear system Ax = b, with A being the /// combined derivative matrix of the residual and b /// being the residual itself. /// \param[in] A matrix A /// \param[inout] x solution to be computed x /// \param[in] b right hand side b template void solve(Operator& opA, Vector& x, Vector& b, Comm& comm) const { Dune::InverseOperatorResult result; // Parallel version is deactivated until we figure out how to do it properly. #if HAVE_MPI if (parallelInformation_.type() == typeid(ParallelISTLInformation)) { const size_t size = opA.getmat().N(); const ParallelISTLInformation& info = boost::any_cast( parallelInformation_); // As we use a dune-istl with block size np the number of components // per parallel is only one. info.copyValuesTo(comm.indexSet(), comm.remoteIndices(), size, 1); // Construct operator, scalar product and vectors needed. constructPreconditionerAndSolve(opA, x, b, comm, result); } else #endif { OPM_THROW(std::logic_error,"this method if for parallel solve only"); } checkConvergence( result ); } /// Solve the linear system Ax = b, with A being the /// combined derivative matrix of the residual and b /// being the residual itself. /// \param[in] A matrix A /// \param[inout] x solution to be computed x /// \param[in] b right hand side b template void solve(Operator& opA, Vector& x, Vector& b ) const { Dune::InverseOperatorResult result; // Construct operator, scalar product and vectors needed. Dune::Amg::SequentialInformation info; constructPreconditionerAndSolve(opA, x, b, info, result); checkConvergence( result ); } void checkConvergence( const Dune::InverseOperatorResult& result ) const { // store number of iterations iterations_ = result.iterations; // Check for failure of linear solver. if (!parameters_.ignoreConvergenceFailure_ && !result.converged) { const std::string msg("Convergence failure for linear solver."); OPM_THROW_NOLOG(LinearSolverProblem, msg); } } protected: mutable int iterations_; boost::any parallelInformation_; bool isIORank_; NewtonIterationBlackoilInterleavedParameters parameters_; }; // end ISTLSolver } // namespace Opm #endif