22#ifndef OPM_PRECONDITIONERFACTORY_HEADER
23#define OPM_PRECONDITIONERFACTORY_HEADER
24#include <opm/common/TimingMacros.hpp>
25#include <opm/simulators/linalg/PreconditionerWithUpdate.hpp>
27#include <dune/istl/paamg/aggregates.hh>
28#include <dune/istl/paamg/matrixhierarchy.hh>
41template <
class Operator,
class Comm,
class Matrix,
class Vector>
44 using PrecPtr = std::shared_ptr<Dune::PreconditionerWithUpdate<Vector, Vector>>;
46 = Dune::Amg::AggregationCriterion<Dune::Amg::SymmetricDependency<Matrix, Dune::Amg::FirstDiagonal>>;
47 using Criterion = Dune::Amg::CoarsenCriterion<CriterionBase>;
51 template <
class Smoother>
52 static PrecPtr makeAmgPreconditioner(
const Operator& op,
54 bool useKamg =
false);
62template <
class Operator,
class Comm>
63class PreconditionerFactory
67 using Matrix =
typename Operator::matrix_type;
68 using Vector =
typename Operator::domain_type;
71 using PrecPtr = std::shared_ptr<Dune::PreconditionerWithUpdate<Vector, Vector>>;
75 const std::function<Vector()>&, std::size_t)>;
77 const std::function<Vector()>&, std::size_t,
const Comm&)>;
86 const std::function<Vector()>& weightsCalculator = {},
87 std::size_t pressureIndex = std::numeric_limits<std::size_t>::max());
97 const std::function<Vector()>& weightsCalculator,
const Comm& comm,
98 std::size_t pressureIndex = std::numeric_limits<std::size_t>::max());
107 std::size_t pressureIndex = std::numeric_limits<std::size_t>::max());
125 static void addCreator(
const std::string& type, ParCreator creator);
128 = Dune::Amg::AggregationCriterion<Dune::Amg::SymmetricDependency<Matrix, Dune::Amg::FirstDiagonal>>;
129 using Criterion = Dune::Amg::CoarsenCriterion<CriterionBase>;
135 static PreconditionerFactory& instance();
138 PreconditionerFactory();
142 const std::function<Vector()> weightsCalculator,
143 std::size_t pressureIndex);
146 const std::function<Vector()> weightsCalculator,
147 std::size_t pressureIndex,
const Comm& comm);
150 void doAddCreator(
const std::string& type,
Creator c);
153 void doAddCreator(
const std::string& type, ParCreator c);
156 std::map<std::string, Creator> creators_;
157 std::map<std::string, ParCreator> parallel_creators_;
158 bool defAdded_=
false;
static void addCreator(const std::string &type, ParCreator creator)
Add a creator for a parallel preconditioner to the PreconditionerFactory.
Definition PreconditionerFactory_impl.hpp:193
std::shared_ptr< Dune::PreconditionerWithUpdate< Vector, Vector > > PrecPtr
Definition PreconditionerFactory.hpp:71
static PrecPtr create(const Operator &op, const PropertyTree &prm, const std::function< Vector()> &weightsCalculator={}, std::size_t pressureIndex=std::numeric_limits< std::size_t >::max())
Create a new serial preconditioner and return a pointer to it.
Definition PreconditionerFactory_impl.hpp:154
static PrecPtr create(const Operator &op, const PropertyTree &prm, const Comm &comm, std::size_t pressureIndex=std::numeric_limits< std::size_t >::max())
Create a new parallel preconditioner and return a pointer to it.
Definition PreconditionerFactory_impl.hpp:176
static void addCreator(const std::string &type, Creator creator)
Add a creator for a serial preconditioner to the PreconditionerFactory.
Definition PreconditionerFactory_impl.hpp:186
std::function< PrecPtr(const OperatorType &, const PropertyTree &, const std::function< Vector()> &, std::size_t)> Creator
Definition PreconditionerFactory.hpp:74
static PrecPtr create(const Operator &op, const PropertyTree &prm, const std::function< Vector()> &weightsCalculator, const Comm &comm, std::size_t pressureIndex=std::numeric_limits< std::size_t >::max())
Create a new parallel preconditioner and return a pointer to it.
Definition PreconditionerFactory_impl.hpp:164
typename OperatorType::matrix_type Matrix
Definition PreconditionerFactory.hpp:67
Hierarchical collection of key/value pairs.
Definition PropertyTree.hpp:39
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45
Definition PreconditionerFactory.hpp:43