/* Copyright 2019 SINTEF Digital, Mathematics and Cybernetics. 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_FLEXIBLE_SOLVER_HEADER_INCLUDED #define OPM_FLEXIBLE_SOLVER_HEADER_INCLUDED #include #include #include #include #include #include namespace Dune { template class FlexibleSolver : Dune::InverseOperator { public: using MatrixType = MatrixTypeT; using VectorType = VectorTypeT; FlexibleSolver(const boost::property_tree::ptree& prm, const MatrixType& matrix) { makeSolver(prm, matrix); } virtual void apply(VectorType& x, VectorType& rhs, Dune::InverseOperatorResult& res) override { linsolver_->apply(x, rhs, res); } virtual void apply(VectorType& x, VectorType& rhs, double reduction, Dune::InverseOperatorResult& res) override { linsolver_->apply(x, rhs, reduction, res); } virtual Dune::SolverCategory::Category category() const override { return Dune::SolverCategory::sequential; } private: void makeSolver(const boost::property_tree::ptree& prm, const MatrixType& matrix) { const double tol = prm.get("tol"); const int maxiter = prm.get("maxiter"); linearoperator_.reset(new Dune::MatrixAdapter(matrix)); preconditioner_ = Dune::makePreconditioner(*linearoperator_, prm); int verbosity = prm.get("verbosity"); std::string solver_type = prm.get("solver"); if (solver_type == "bicgstab") { linsolver_.reset(new Dune::BiCGSTABSolver(*linearoperator_, *preconditioner_, tol, // desired residual reduction factor maxiter, // maximum number of iterations verbosity)); } else if (solver_type == "loopsolver") { linsolver_.reset(new Dune::LoopSolver(*linearoperator_, *preconditioner_, tol, // desired residual reduction factor maxiter, // maximum number of iterations verbosity)); } else if (solver_type == "gmres") { int restart = prm.get("restart"); linsolver_.reset(new Dune::RestartedGMResSolver(*linearoperator_, *preconditioner_, tol, restart, // desired residual reduction factor maxiter, // maximum number of iterations verbosity)); #if HAVE_SUITESPARSE_UMFPACK } else if (solver_type == "umfpack") { bool dummy = false; linsolver_.reset(new Dune::UMFPack(linearoperator_->getmat(), verbosity, dummy)); #endif } else { std::string msg("Solver not known "); msg += solver_type; throw std::runtime_error(msg); } } std::shared_ptr> preconditioner_; std::shared_ptr> linearoperator_; std::shared_ptr> linsolver_; }; } // namespace Dune #endif // OPM_FLEXIBLE_SOLVER_HEADER_INCLUDED