2019-05-20 06:23:57 -05:00
|
|
|
/*
|
|
|
|
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 <http://www.gnu.org/licenses/>.
|
|
|
|
*/
|
|
|
|
|
|
|
|
#ifndef OPM_OWNINGTWOLEVELPRECONDITIONER_HEADER_INCLUDED
|
|
|
|
#define OPM_OWNINGTWOLEVELPRECONDITIONER_HEADER_INCLUDED
|
|
|
|
|
2019-05-20 07:31:36 -05:00
|
|
|
#include <opm/simulators/linalg/PreconditionerWithUpdate.hpp>
|
2019-05-20 06:23:57 -05:00
|
|
|
#include <opm/simulators/linalg/PressureSolverPolicy.hpp>
|
|
|
|
#include <opm/simulators/linalg/PressureTransferPolicy.hpp>
|
|
|
|
#include <opm/simulators/linalg/getQuasiImpesWeights.hpp>
|
2019-05-20 07:31:36 -05:00
|
|
|
#include <opm/simulators/linalg/twolevelmethodcpr.hh>
|
2019-05-20 06:23:57 -05:00
|
|
|
|
|
|
|
#include <dune/common/fmatrix.hh>
|
|
|
|
#include <dune/istl/bcrsmatrix.hh>
|
|
|
|
#include <dune/istl/paamg/amg.hh>
|
|
|
|
|
|
|
|
#include <boost/property_tree/ptree.hpp>
|
|
|
|
|
|
|
|
#include <fstream>
|
2019-05-20 07:31:36 -05:00
|
|
|
#include <type_traits>
|
2019-05-20 06:23:57 -05:00
|
|
|
|
|
|
|
|
2020-01-07 10:53:54 -06:00
|
|
|
namespace Opm
|
|
|
|
{
|
2019-05-29 09:21:34 -05:00
|
|
|
// Circular dependency between PreconditionerFactory [which can make an OwningTwoLevelPreconditioner]
|
|
|
|
// and OwningTwoLevelPreconditioner [which uses PreconditionerFactory to choose the fine-level smoother]
|
2019-05-20 06:23:57 -05:00
|
|
|
// must be broken, accomplished by forward-declaration here.
|
2019-05-29 09:21:34 -05:00
|
|
|
template <class Operator, class Comm = Dune::Amg::SequentialInformation>
|
|
|
|
class PreconditionerFactory;
|
2020-01-07 10:53:54 -06:00
|
|
|
}
|
|
|
|
|
|
|
|
namespace Dune
|
|
|
|
{
|
2019-05-20 07:31:36 -05:00
|
|
|
|
2019-05-20 06:23:57 -05:00
|
|
|
|
|
|
|
// Must forward-declare FlexibleSolver as we want to use it as solver for the pressure system.
|
|
|
|
template <class MatrixTypeT, class VectorTypeT>
|
|
|
|
class FlexibleSolver;
|
|
|
|
|
2019-05-29 09:21:34 -05:00
|
|
|
|
|
|
|
/// A version of the two-level preconditioner that is:
|
|
|
|
/// - Self-contained, because it owns its policy components.
|
|
|
|
/// - Flexible, because it uses the runtime-flexible solver
|
|
|
|
/// and preconditioner factory.
|
2019-05-20 07:31:36 -05:00
|
|
|
template <class OperatorType,
|
|
|
|
class VectorType,
|
|
|
|
bool transpose = false,
|
|
|
|
class Communication = Dune::Amg::SequentialInformation>
|
|
|
|
class OwningTwoLevelPreconditioner : public Dune::PreconditionerWithUpdate<VectorType, VectorType>
|
2019-05-20 06:23:57 -05:00
|
|
|
{
|
|
|
|
public:
|
|
|
|
using pt = boost::property_tree::ptree;
|
2019-05-20 07:31:36 -05:00
|
|
|
using MatrixType = typename OperatorType::matrix_type;
|
2020-01-07 10:53:54 -06:00
|
|
|
using PrecFactory = Opm::PreconditionerFactory<OperatorType, Communication>;
|
2019-05-20 06:23:57 -05:00
|
|
|
|
2019-05-28 09:22:54 -05:00
|
|
|
OwningTwoLevelPreconditioner(const OperatorType& linearoperator, const pt& prm)
|
2019-05-20 07:31:36 -05:00
|
|
|
: linear_operator_(linearoperator)
|
2019-05-29 09:21:34 -05:00
|
|
|
, finesmoother_(PrecFactory::create(linearoperator, prm.get_child("finesmoother")))
|
2019-05-20 07:31:36 -05:00
|
|
|
, comm_(nullptr)
|
2019-05-20 06:23:57 -05:00
|
|
|
, weights_(Opm::Amg::getQuasiImpesWeights<MatrixType, VectorType>(
|
|
|
|
linearoperator.getmat(), prm.get<int>("pressure_var_index"), transpose))
|
2019-08-07 07:17:17 -05:00
|
|
|
, levelTransferPolicy_(dummy_comm_, weights_, prm.get<int>("pressure_var_index"))
|
2019-05-20 06:23:57 -05:00
|
|
|
, coarseSolverPolicy_(prm.get_child("coarsesolver"))
|
|
|
|
, twolevel_method_(linearoperator,
|
|
|
|
finesmoother_,
|
|
|
|
levelTransferPolicy_,
|
|
|
|
coarseSolverPolicy_,
|
|
|
|
transpose ? 1 : 0,
|
|
|
|
transpose ? 0 : 1)
|
2019-05-20 07:31:36 -05:00
|
|
|
, prm_(prm)
|
|
|
|
{
|
|
|
|
if (prm.get<int>("verbosity") > 10) {
|
|
|
|
std::ofstream outfile(prm.get<std::string>("weights_filename"));
|
|
|
|
if (!outfile) {
|
|
|
|
throw std::runtime_error("Could not write weights");
|
|
|
|
}
|
|
|
|
Dune::writeMatrixMarket(weights_, outfile);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2019-05-28 09:22:54 -05:00
|
|
|
OwningTwoLevelPreconditioner(const OperatorType& linearoperator, const pt& prm, const Communication& comm)
|
2019-05-20 07:31:36 -05:00
|
|
|
: linear_operator_(linearoperator)
|
2019-05-29 09:21:34 -05:00
|
|
|
, finesmoother_(PrecFactory::create(linearoperator, prm.get_child("finesmoother"), comm))
|
2019-05-20 07:31:36 -05:00
|
|
|
, comm_(&comm)
|
|
|
|
, weights_(Opm::Amg::getQuasiImpesWeights<MatrixType, VectorType>(
|
|
|
|
linearoperator.getmat(), prm.get<int>("pressure_var_index"), transpose))
|
|
|
|
, levelTransferPolicy_(*comm_, weights_, prm.get<int>("pressure_var_index"))
|
|
|
|
, coarseSolverPolicy_(prm.get_child("coarsesolver"))
|
|
|
|
, twolevel_method_(linearoperator,
|
|
|
|
finesmoother_,
|
|
|
|
levelTransferPolicy_,
|
|
|
|
coarseSolverPolicy_,
|
|
|
|
transpose ? 1 : 0,
|
|
|
|
transpose ? 0 : 1)
|
|
|
|
, prm_(prm)
|
2019-05-20 06:23:57 -05:00
|
|
|
{
|
|
|
|
if (prm.get<int>("verbosity") > 10) {
|
|
|
|
std::ofstream outfile(prm.get<std::string>("weights_filename"));
|
|
|
|
if (!outfile) {
|
|
|
|
throw std::runtime_error("Could not write weights");
|
|
|
|
}
|
|
|
|
Dune::writeMatrixMarket(weights_, outfile);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
virtual void pre(VectorType& x, VectorType& b) override
|
|
|
|
{
|
|
|
|
twolevel_method_.pre(x, b);
|
|
|
|
}
|
2019-05-20 07:31:36 -05:00
|
|
|
|
2019-05-20 06:23:57 -05:00
|
|
|
virtual void apply(VectorType& v, const VectorType& d) override
|
|
|
|
{
|
|
|
|
twolevel_method_.apply(v, d);
|
|
|
|
}
|
2019-05-20 07:31:36 -05:00
|
|
|
|
2019-05-20 06:23:57 -05:00
|
|
|
virtual void post(VectorType& x) override
|
|
|
|
{
|
|
|
|
twolevel_method_.post(x);
|
|
|
|
}
|
2019-05-20 07:31:36 -05:00
|
|
|
|
|
|
|
virtual void update() override
|
|
|
|
{
|
|
|
|
Opm::Amg::getQuasiImpesWeights<MatrixType, VectorType>(
|
|
|
|
linear_operator_.getmat(), prm_.get<int>("pressure_var_index"), transpose, weights_);
|
2019-06-06 03:36:17 -05:00
|
|
|
updateImpl(comm_);
|
2019-05-20 07:31:36 -05:00
|
|
|
}
|
|
|
|
|
2019-05-20 06:23:57 -05:00
|
|
|
virtual Dune::SolverCategory::Category category() const override
|
|
|
|
{
|
2019-05-20 07:31:36 -05:00
|
|
|
return linear_operator_.category();
|
2019-05-20 06:23:57 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
private:
|
|
|
|
using PressureMatrixType = Dune::BCRSMatrix<Dune::FieldMatrix<double, 1, 1>>;
|
|
|
|
using PressureVectorType = Dune::BlockVector<Dune::FieldVector<double, 1>>;
|
2019-05-20 07:31:36 -05:00
|
|
|
using SeqCoarseOperatorType = Dune::MatrixAdapter<PressureMatrixType, PressureVectorType, PressureVectorType>;
|
|
|
|
using ParCoarseOperatorType
|
|
|
|
= Dune::OverlappingSchwarzOperator<PressureMatrixType, PressureVectorType, PressureVectorType, Communication>;
|
|
|
|
using CoarseOperatorType = std::conditional_t<std::is_same<Communication, Dune::Amg::SequentialInformation>::value,
|
|
|
|
SeqCoarseOperatorType,
|
|
|
|
ParCoarseOperatorType>;
|
2019-05-20 06:23:57 -05:00
|
|
|
using LevelTransferPolicy = Opm::PressureTransferPolicy<OperatorType, CoarseOperatorType, Communication, transpose>;
|
2019-05-20 07:31:36 -05:00
|
|
|
using CoarseSolverPolicy = Dune::Amg::PressureSolverPolicy<CoarseOperatorType,
|
|
|
|
FlexibleSolver<PressureMatrixType, PressureVectorType>,
|
|
|
|
LevelTransferPolicy>;
|
2019-05-20 06:23:57 -05:00
|
|
|
using TwoLevelMethod
|
2019-05-20 07:31:36 -05:00
|
|
|
= Dune::Amg::TwoLevelMethodCpr<OperatorType, CoarseSolverPolicy, Dune::Preconditioner<VectorType, VectorType>>;
|
|
|
|
|
2019-05-29 09:21:34 -05:00
|
|
|
// Handling parallel vs serial instantiation of preconditioner factory.
|
2019-05-20 07:31:36 -05:00
|
|
|
template <class Comm>
|
2019-06-06 03:36:17 -05:00
|
|
|
void updateImpl(const Comm*)
|
2019-05-20 07:31:36 -05:00
|
|
|
{
|
|
|
|
// Parallel case.
|
2019-05-29 09:21:34 -05:00
|
|
|
finesmoother_ = PrecFactory::create(linear_operator_, prm_.get_child("finesmoother"), *comm_);
|
2019-05-20 07:31:36 -05:00
|
|
|
twolevel_method_.updatePreconditioner(finesmoother_, coarseSolverPolicy_);
|
|
|
|
}
|
|
|
|
|
2019-06-06 03:36:17 -05:00
|
|
|
void updateImpl(const Dune::Amg::SequentialInformation*)
|
2019-05-20 07:31:36 -05:00
|
|
|
{
|
|
|
|
// Serial case.
|
2019-05-29 09:21:34 -05:00
|
|
|
finesmoother_ = PrecFactory::create(linear_operator_, prm_.get_child("finesmoother"));
|
2019-05-20 07:31:36 -05:00
|
|
|
twolevel_method_.updatePreconditioner(finesmoother_, coarseSolverPolicy_);
|
|
|
|
}
|
2019-05-20 06:23:57 -05:00
|
|
|
|
2019-05-28 09:22:54 -05:00
|
|
|
const OperatorType& linear_operator_;
|
2019-05-20 06:23:57 -05:00
|
|
|
std::shared_ptr<Dune::Preconditioner<VectorType, VectorType>> finesmoother_;
|
2019-05-20 07:31:36 -05:00
|
|
|
const Communication* comm_;
|
2019-05-20 06:23:57 -05:00
|
|
|
VectorType weights_;
|
|
|
|
LevelTransferPolicy levelTransferPolicy_;
|
|
|
|
CoarseSolverPolicy coarseSolverPolicy_;
|
|
|
|
TwoLevelMethod twolevel_method_;
|
2019-05-20 07:31:36 -05:00
|
|
|
boost::property_tree::ptree prm_;
|
2019-08-07 07:17:17 -05:00
|
|
|
Communication dummy_comm_;
|
2019-05-20 06:23:57 -05:00
|
|
|
};
|
|
|
|
|
|
|
|
} // namespace Dune
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
#endif // OPM_OWNINGTWOLEVELPRECONDITIONER_HEADER_INCLUDED
|