1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177
|
// -*- 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 <http://www.gnu.org/licenses/>.
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::ParallelBiCGStabSolverBackend
*/
#ifndef EWOMS_PARALLEL_BICGSTAB_BACKEND_HH
#define EWOMS_PARALLEL_BICGSTAB_BACKEND_HH
#include "linalgproperties.hh"
#include "parallelbasebackend.hh"
#include "bicgstabsolver.hh"
#include "combinedcriterion.hh"
#include "istlsparsematrixadapter.hh"
#include <memory>
namespace Opm::Linear {
template <class TypeTag>
class ParallelBiCGStabSolverBackend;
} // namespace Opm::Linear
namespace Opm::Properties {
// Create new type tags
namespace TTag {
struct ParallelBiCGStabLinearSolver { using InheritsFrom = std::tuple<ParallelBaseLinearSolver>; };
} // end namespace TTag
template<class TypeTag>
struct LinearSolverBackend<TypeTag, TTag::ParallelBiCGStabLinearSolver>
{ using type = Opm::Linear::ParallelBiCGStabSolverBackend<TypeTag>; };
template<class TypeTag>
struct LinearSolverMaxError<TypeTag, TTag::ParallelBiCGStabLinearSolver>
{
using type = GetPropType<TypeTag, Scalar>;
static constexpr type value = 1e7;
};
} // namespace Opm::Properties
namespace Opm {
namespace Linear {
/*!
* \ingroup Linear
*
* \brief Implements a generic linear solver abstraction.
*
* Chosing the preconditioner works by setting the "PreconditionerWrapper" property:
*
* \code
* template<class TypeTag>
* struct PreconditionerWrapper<TypeTag, TTag::YourTypeTag>
* { using type = Opm::Linear::PreconditionerWrapper$PRECONDITIONER<TypeTag>; };
* \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: An ILU(0) preconditioner. The results of this
* preconditioner are the same as setting the
* PreconditionerOrder property to 0 and using the ILU(n)
* preconditioner. The reason for the existence of ILU0 is
* that it is computationally cheaper because it does not
* need to consider things which are only required for
* higher orders
*/
template <class TypeTag>
class ParallelBiCGStabSolverBackend : public ParallelBaseBackend<TypeTag>
{
using ParentType = ParallelBaseBackend<TypeTag>;
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
using Simulator = GetPropType<TypeTag, Properties::Simulator>;
using SparseMatrixAdapter = GetPropType<TypeTag, Properties::SparseMatrixAdapter>;
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 = BiCGStabSolver<ParallelOperator,
OverlappingVector,
ParallelPreconditioner>;
static_assert(std::is_same<SparseMatrixAdapter, IstlSparseMatrixAdapter<MatrixBlock> >::value,
"The ParallelIstlSolverBackend linear solver backend requires the IstlSparseMatrixAdapter");
public:
ParallelBiCGStabSolverBackend(const Simulator& simulator)
: ParentType(simulator)
{ }
static void registerParameters()
{
ParentType::registerParameters();
EWOMS_REGISTER_PARAM(TypeTag, Scalar, LinearSolverMaxError,
"The maximum residual error which the linear solver tolerates"
" without giving up");
}
protected:
friend ParentType;
std::shared_ptr<RawLinearSolver> prepareSolver_(ParallelOperator& parOperator,
ParallelScalarProduct& parScalarProduct,
ParallelPreconditioner& parPreCond)
{
const auto& gridView = this->simulator_.gridView();
using CCC = CombinedCriterion<OverlappingVector, decltype(gridView.comm())>;
Scalar linearSolverTolerance = EWOMS_GET_PARAM(TypeTag, Scalar, LinearSolverTolerance);
Scalar linearSolverAbsTolerance = EWOMS_GET_PARAM(TypeTag, Scalar, LinearSolverAbsTolerance);
if(linearSolverAbsTolerance < 0.0)
linearSolverAbsTolerance = this->simulator_.model().newtonMethod().tolerance() / 100.0;
convCrit_.reset(new CCC(gridView.comm(),
/*residualReductionTolerance=*/linearSolverTolerance,
/*absoluteResidualTolerance=*/linearSolverAbsTolerance,
EWOMS_GET_PARAM(TypeTag, Scalar, LinearSolverMaxError)));
auto bicgstabSolver =
std::make_shared<RawLinearSolver>(parPreCond, *convCrit_, parScalarProduct);
int verbosity = 0;
if (parOperator.overlap().myRank() == 0)
verbosity = EWOMS_GET_PARAM(TypeTag, int, LinearSolverVerbosity);
bicgstabSolver->setVerbosity(verbosity);
bicgstabSolver->setMaxIterations(EWOMS_GET_PARAM(TypeTag, int, LinearSolverMaxIterations));
bicgstabSolver->setLinearOperator(&parOperator);
bicgstabSolver->setRhs(this->overlappingb_);
return bicgstabSolver;
}
std::pair<bool,int> runSolver_(std::shared_ptr<RawLinearSolver> solver)
{
bool converged = solver->apply(*this->overlappingx_);
return std::make_pair(converged, int(solver->report().iterations()));
}
void cleanupSolver_()
{ /* nothing to do */ }
std::unique_ptr<ConvergenceCriterion<OverlappingVector> > convCrit_;
};
}} // namespace Linear, Opm
#endif
|