/*
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
#include
#include
namespace Dune
{
template
struct IsComm : std::false_type
{};
template<>
struct IsComm : std::true_type
{};
#if HAVE_MPI
template
struct IsComm> : std::true_type
{};
#endif
/// A solver class that encapsulates all needed objects for a linear solver
/// (operator, scalar product, iterative solver and preconditioner) and sets
/// them up based on runtime parameters, using the PreconditionerFactory for
/// setting up preconditioners.
template
class FlexibleSolver : public Dune::InverseOperator
{
public:
using MatrixType = MatrixTypeT;
using VectorType = VectorTypeT;
/// Create a sequential solver.
FlexibleSolver(const boost::property_tree::ptree& prm, const MatrixType& matrix,
const std::function& weightsCalculator = std::function())
{
init(prm, matrix, weightsCalculator, Dune::Amg::SequentialInformation());
}
/// Create a parallel solver (if Comm is e.g. OwnerOverlapCommunication).
template
FlexibleSolver(const boost::property_tree::ptree& prm,
const MatrixType& matrix,
const typename std::enable_if::value, Comm>::type& comm)
{
init(prm, matrix, std::function(), comm);
}
/// Create a parallel solver (if Comm is e.g. OwnerOverlapCommunication).
template
FlexibleSolver(const boost::property_tree::ptree& prm, const MatrixType& matrix,
const std::function& weightsCalculator, const Comm& comm)
{
init(prm, matrix, weightsCalculator, comm);
}
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);
}
/// Type of the contained preconditioner.
using AbstractPrecondType = Dune::PreconditionerWithUpdate;
/// Access the contained preconditioner.
AbstractPrecondType& preconditioner()
{
return *preconditioner_;
}
virtual Dune::SolverCategory::Category category() const override
{
return linearoperator_->category();
}
private:
using AbstractOperatorType = Dune::AssembledLinearOperator;
using AbstractScalarProductType = Dune::ScalarProduct;
using AbstractSolverType = Dune::InverseOperator;
// Machinery for making sequential or parallel operators/preconditioners/scalar products.
template
void initOpPrecSp(const MatrixType& matrix, const boost::property_tree::ptree& prm,
const std::function weightsCalculator, const Comm& comm)
{
// Parallel case.
using ParOperatorType = Dune::OverlappingSchwarzOperator;
using pt = const boost::property_tree::ptree;
auto linop = std::make_shared(matrix, comm);
linearoperator_ = linop;
auto child = prm.get_child_optional("preconditioner");
preconditioner_
= Opm::PreconditionerFactory::create(*linop, child? *child : pt(),
weightsCalculator, comm);
scalarproduct_ = Dune::createScalarProduct(comm, linearoperator_->category());
}
void initOpPrecSp(const MatrixType& matrix, const boost::property_tree::ptree& prm,
const std::function weightsCalculator, const Dune::Amg::SequentialInformation&)
{
// Sequential case.
using SeqOperatorType = Dune::MatrixAdapter;
using pt = const boost::property_tree::ptree;
auto linop = std::make_shared(matrix);
linearoperator_ = linop;
auto child = prm.get_child_optional("preconditioner");
preconditioner_ = Opm::PreconditionerFactory::create(*linop, child? *child : pt(),
weightsCalculator);
scalarproduct_ = std::make_shared>();
}
void initSolver(const boost::property_tree::ptree& prm)
{
const double tol = prm.get("tol", 1e-2);
const int maxiter = prm.get("maxiter", 200);
const int verbosity = prm.get("verbosity", 0);
const std::string solver_type = prm.get("solver", "bicgstab");
if (solver_type == "bicgstab") {
linsolver_.reset(new Dune::BiCGSTABSolver(*linearoperator_,
*scalarproduct_,
*preconditioner_,
tol, // desired residual reduction factor
maxiter, // maximum number of iterations
verbosity));
} else if (solver_type == "loopsolver") {
linsolver_.reset(new Dune::LoopSolver(*linearoperator_,
*scalarproduct_,
*preconditioner_,
tol, // desired residual reduction factor
maxiter, // maximum number of iterations
verbosity));
} else if (solver_type == "gmres") {
int restart = prm.get("restart", 15);
linsolver_.reset(new Dune::RestartedGMResSolver(*linearoperator_,
*scalarproduct_,
*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 {
OPM_THROW(std::invalid_argument, "Properties: Solver " << solver_type << " not known.");
}
}
// Main initialization routine.
// Call with Comm == Dune::Amg::SequentialInformation to get a serial solver.
template
void init(const boost::property_tree::ptree& prm, const MatrixType& matrix,
const std::function weightsCalculator, const Comm& comm)
{
initOpPrecSp(matrix, prm, weightsCalculator, comm);
initSolver(prm);
}
std::shared_ptr linearoperator_;
std::shared_ptr preconditioner_;
std::shared_ptr scalarproduct_;
std::shared_ptr linsolver_;
};
} // namespace Dune
#endif // OPM_FLEXIBLE_SOLVER_HEADER_INCLUDED