// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*- // vi: set et ts=4 sw=4 sts=4: /* 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 2 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 . Consult the COPYING file in the top-level source directory of this module for the precise wording of the license and the list of copyright holders. */ /*! * \file * \copydoc Opm::Linear::ParallelIstlSolverBackend */ #ifndef EWOMS_PARALLEL_ISTL_BACKEND_HH #define EWOMS_PARALLEL_ISTL_BACKEND_HH #include "linalgproperties.hh" #include "parallelbasebackend.hh" #include "istlsolverwrappers.hh" #include "istlsparsematrixadapter.hh" #include namespace Opm::Properties::TTag { // Create new type tag struct ParallelIstlLinearSolver { using InheritsFrom = std::tuple; }; } // namespace Opm::Properties::TTag namespace Opm::Linear { /*! * \ingroup Linear * * \brief Provides all unmodified linear solvers from dune-istl * * To set the linear solver, use * \code * template * struct LinearSolverWrapper * { using type = Opm::Linear::SolverWrapper$SOLVER; }; * \endcode * * The possible choices for '\c $SOLVER' are: * - \c Richardson: A fixpoint solver using the Richardson iteration * - \c SteepestDescent: The steepest descent solver * - \c ConjugatedGradients: A conjugated gradients solver * - \c BiCGStab: A stabilized bi-conjugated gradients solver * - \c MinRes: A solver based on the minimized residual algorithm * - \c RestartedGMRes: A restarted GMRES solver * * Chosing the preconditioner works in an analogous way: * \code * template * struct PreconditionerWrapper * { using type = Opm::Linear::PreconditionerWrapper$PRECONDITIONER; }; * \endcode * * Where the choices possible for '\c $PRECONDITIONER' are: * - \c Jacobi: A Jacobi preconditioner * - \c GaussSeidel: A Gauss-Seidel preconditioner * - \c SSOR: A symmetric successive overrelaxation (SSOR) preconditioner * - \c SOR: A successive overrelaxation (SOR) preconditioner * - \c ILUn: An ILU(n) preconditioner * - \c ILU0: A specialized (and optimized) ILU(0) preconditioner */ template class ParallelIstlSolverBackend : public ParallelBaseBackend { using ParentType = ParallelBaseBackend; using Scalar = GetPropType; using Simulator = GetPropType; using LinearSolverWrapper = GetPropType; using SparseMatrixAdapter = GetPropType; using ParallelOperator = typename ParentType::ParallelOperator; using OverlappingVector = typename ParentType::OverlappingVector; using ParallelPreconditioner = typename ParentType::ParallelPreconditioner; using ParallelScalarProduct = typename ParentType::ParallelScalarProduct; using MatrixBlock = typename SparseMatrixAdapter::MatrixBlock; using RawLinearSolver = typename LinearSolverWrapper::RawSolver; static_assert(std::is_same >::value, "The ParallelIstlSolverBackend linear solver backend requires the IstlSparseMatrixAdapter"); public: ParallelIstlSolverBackend(const Simulator& simulator) : ParentType(simulator) { } /*! * \brief Register all run-time parameters for the linear solver. */ static void registerParameters() { ParentType::registerParameters(); LinearSolverWrapper::registerParameters(); } protected: friend ParentType; std::shared_ptr prepareSolver_(ParallelOperator& parOperator, ParallelScalarProduct& parScalarProduct, ParallelPreconditioner& parPreCond) { return solverWrapper_.get(parOperator, parScalarProduct, parPreCond); } void cleanupSolver_() { solverWrapper_.cleanup(); } std::pair runSolver_(std::shared_ptr solver) { Dune::InverseOperatorResult result; solver->apply(*this->overlappingx_, *this->overlappingb_, result); return std::make_pair(result.converged, result.iterations); } LinearSolverWrapper solverWrapper_; }; } // namespace Opm::Linear namespace Opm::Properties { template struct LinearSolverBackend { using type = Opm::Linear::ParallelIstlSolverBackend; }; template struct LinearSolverWrapper { using type = Opm::Linear::SolverWrapperBiCGStab; }; template struct PreconditionerWrapper { using type = Opm::Linear::PreconditionerWrapperILU; }; //! set the GMRes restart parameter to 10 by default template struct GMResRestart { static constexpr int value = 10; }; } // namespace Opm::Properties #endif