PreconditionerFactory_impl.hpp
Go to the documentation of this file.
1/*
2 Copyright 2009, 2010 SINTEF ICT, Applied Mathematics.
3 Copyright 2019 SINTEF Digital, Mathematics and Cybernetics.
4
5 This file is part of the Open Porous Media project (OPM).
6
7 OPM is free software: you can redistribute it and/or modify
8 it under the terms of the GNU General Public License as published by
9 the Free Software Foundation, either version 3 of the License, or
10 (at your option) any later version.
11
12 OPM is distributed in the hope that it will be useful,
13 but WITHOUT ANY WARRANTY; without even the implied warranty of
14 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 GNU General Public License for more details.
16
17 You should have received a copy of the GNU General Public License
18 along with OPM. If not, see <http://www.gnu.org/licenses/>.
19*/
20
21#include <config.h>
22
23#include <opm/common/ErrorMacros.hpp>
24#include <opm/common/TimingMacros.hpp>
25
27
42
43#include <dune/common/unused.hh>
44#include <dune/istl/owneroverlapcopy.hh>
45#include <dune/istl/paamg/amg.hh>
46#include <dune/istl/paamg/fastamg.hh>
47#include <dune/istl/paamg/kamg.hh>
48#include <dune/istl/preconditioners.hh>
49
50// Include all cuistl/GPU preconditioners inside of this headerfile
52
53#if HAVE_HYPRE
55#endif
56
57#if HAVE_AMGX
59#endif
60
61#include <cassert>
62
64
65namespace Opm {
66
67
68
69template <class Operator, class Comm>
70PreconditionerFactory<Operator, Comm>::PreconditionerFactory()
71{
72}
73
74
75template <class Operator, class Comm>
76PreconditionerFactory<Operator, Comm>&
77PreconditionerFactory<Operator, Comm>::instance()
78{
79 static PreconditionerFactory singleton;
80 return singleton;
81}
82
83template <class Operator, class Comm>
85PreconditionerFactory<Operator, Comm>::doCreate(const Operator& op,
86 const PropertyTree& prm,
87 const std::function<Vector()> weightsCalculator,
88 std::size_t pressureIndex)
89{
90 if (!defAdded_) {
92 defAdded_ = true;
93 }
94 const std::string& type = prm.get<std::string>("type", "ParOverILU0");
95 auto it = creators_.find(type);
96 if (it == creators_.end()) {
97 std::ostringstream msg;
98 msg << "Preconditioner type " << type << " is not registered in the factory. Available types are: ";
99 for (const auto& prec : creators_) {
100 msg << prec.first << ' ';
101 }
102 msg << std::endl;
103 OPM_THROW(std::invalid_argument, msg.str());
104 }
105 return it->second(op, prm, weightsCalculator, pressureIndex);
106}
107
108template <class Operator, class Comm>
110PreconditionerFactory<Operator, Comm>::doCreate(const Operator& op,
111 const PropertyTree& prm,
112 const std::function<Vector()> weightsCalculator,
113 std::size_t pressureIndex,
114 const Comm& comm)
115{
116 if (!defAdded_) {
118 defAdded_ = true;
119 }
120 const std::string& type = prm.get<std::string>("type", "ParOverILU0");
121 auto it = parallel_creators_.find(type);
122 if (it == parallel_creators_.end()) {
123 std::ostringstream msg;
124 msg << "Parallel preconditioner type " << type << " is not registered in the factory. Available types are: ";
125 for (const auto& prec : parallel_creators_) {
126 msg << prec.first << ' ';
127 }
128 msg << std::endl;
129 OPM_THROW(std::invalid_argument, msg.str());
130 }
131 return it->second(op, prm, weightsCalculator, pressureIndex, comm);
132}
133
134template <class Operator, class Comm>
135void
136PreconditionerFactory<Operator, Comm>::doAddCreator(const std::string& type, Creator c)
137{
138 creators_[type] = c;
139}
140
141template <class Operator, class Comm>
142void
143PreconditionerFactory<Operator, Comm>::doAddCreator(const std::string& type, ParCreator c)
144{
145 parallel_creators_[type] = c;
146}
147
148template <class Operator, class Comm>
151 const PropertyTree& prm,
152 const std::function<Vector()>& weightsCalculator,
153 std::size_t pressureIndex)
154{
155 return instance().doCreate(op, prm, weightsCalculator, pressureIndex);
156}
157
158template <class Operator, class Comm>
161 const PropertyTree& prm,
162 const std::function<Vector()>& weightsCalculator,
163 const Comm& comm,
164 std::size_t pressureIndex)
165{
166 return instance().doCreate(op, prm, weightsCalculator, pressureIndex, comm);
167}
168
169
170template <class Operator, class Comm>
173 const PropertyTree& prm,
174 const Comm& comm,
175 std::size_t pressureIndex)
176{
177 return instance().doCreate(op, prm, std::function<Vector()>(), pressureIndex, comm);
178}
179
180template <class Operator, class Comm>
181void
183{
184 instance().doAddCreator(type, creator);
185}
186
187template <class Operator, class Comm>
188void
190{
191 instance().doAddCreator(type, creator);
192}
193
194using CommSeq = Dune::Amg::SequentialInformation;
195
196template<class Scalar, int Dim>
197using OpFSeq = Dune::MatrixAdapter<Dune::BCRSMatrix<Dune::FieldMatrix<Scalar, Dim, Dim>>,
198 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>,
199 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>>;
200template<class Scalar, int Dim>
201using OpBSeq = Dune::MatrixAdapter<Dune::BCRSMatrix<MatrixBlock<Scalar, Dim, Dim>>,
202 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>,
203 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>>;
204
205template<class Scalar, int Dim, bool overlap>
207 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>,
208 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>>;
209
210template<class Scalar, int Dim, bool overlap>
212 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>,
213 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>,
214 overlap>;
215
216#if HAVE_MPI
217using CommPar = Dune::OwnerOverlapCopyCommunication<int, int>;
218
219template<class Scalar, int Dim>
220using OpFPar = Dune::OverlappingSchwarzOperator<Dune::BCRSMatrix<Dune::FieldMatrix<Scalar, Dim, Dim>>,
221 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>,
222 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>,
223 CommPar>;
224
225template<class Scalar, int Dim>
226using OpBPar = Dune::OverlappingSchwarzOperator<Dune::BCRSMatrix<MatrixBlock<Scalar, Dim, Dim>>,
227 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>,
228 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>,
229 CommPar>;
230template<class Scalar, int Dim>
232 Dune::BlockVector<Dune::FieldVector<Scalar,Dim>>,
233 Dune::BlockVector<Dune::FieldVector<Scalar,Dim>>,
234 CommPar>;
235
236template<class Scalar, int Dim>
238 Dune::BlockVector<Dune::FieldVector<Scalar,Dim>>,
239 Dune::BlockVector<Dune::FieldVector<Scalar,Dim>>,
240 CommPar>;
241
242#define INSTANTIATE_PF_PAR(T,Dim) \
243 template class PreconditionerFactory<OpBSeq<T,Dim>, CommPar>; \
244 template class PreconditionerFactory<OpFPar<T,Dim>, CommPar>; \
245 template class PreconditionerFactory<OpBPar<T,Dim>, CommPar>; \
246 template class PreconditionerFactory<OpGLFPar<T,Dim>, CommPar>; \
247 template class PreconditionerFactory<OpGLBPar<T,Dim>, CommPar>; \
248 template class PreconditionerFactory<OpW<T,Dim, false>, CommPar>; \
249 template class PreconditionerFactory<OpWG<T,Dim, true>, CommPar>; \
250 template class PreconditionerFactory<OpBPar<T,Dim>, CommSeq>; \
251 template class PreconditionerFactory<OpGLBPar<T,Dim>, CommSeq>;
252#endif
253
254#define INSTANTIATE_PF_SEQ(T,Dim) \
255 template class PreconditionerFactory<OpFSeq<T,Dim>, CommSeq>; \
256 template class PreconditionerFactory<OpBSeq<T,Dim>, CommSeq>; \
257 template class PreconditionerFactory<OpW<T,Dim, false>, CommSeq>; \
258 template class PreconditionerFactory<OpWG<T,Dim, true>, CommSeq>;
259
260#if HAVE_MPI
261#define INSTANTIATE_PF(T,Dim) \
262 INSTANTIATE_PF_PAR(T,Dim) \
263 INSTANTIATE_PF_SEQ(T,Dim)
264#else
265#define INSTANTIATE_PF(T,Dim) INSTANTIATE_PF_SEQ(T,Dim)
266#endif
267
268} // namespace Opm
Dune::OwnerOverlapCopyCommunication< int, int > Comm
Definition: FlexibleSolver_impl.hpp:304
The AMG preconditioner.
Dune linear operator that assumes ghost rows are ordered after interior rows. Avoids some computation...
Definition: WellOperators.hpp:402
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())
Definition: PreconditionerFactory_impl.hpp:150
std::function< PrecPtr(const Operator &, const PropertyTree &, const std::function< Vector()> &, std::size_t, const Comm &)> ParCreator
Definition: PreconditionerFactory.hpp:77
typename Operator::domain_type Vector
Definition: PreconditionerFactory.hpp:68
std::function< PrecPtr(const Operator &, const PropertyTree &, const std::function< Vector()> &, std::size_t)> Creator
The type of creator functions passed to addCreator().
Definition: PreconditionerFactory.hpp:75
static void addCreator(const std::string &type, Creator creator)
Definition: PreconditionerFactory_impl.hpp:182
std::shared_ptr< Dune::PreconditionerWithUpdate< Vector, Vector > > PrecPtr
The type of pointer returned by create().
Definition: PreconditionerFactory.hpp:71
Hierarchical collection of key/value pairs.
Definition: PropertyTree.hpp:39
Adapter to combine a matrix and another linear operator into a combined linear operator.
Definition: WellOperators.hpp:299
Adapter to combine a matrix and another linear operator into a combined linear operator.
Definition: WellOperators.hpp:225
Definition: blackoilboundaryratevector.hh:39
Dune::MatrixAdapter< Dune::BCRSMatrix< MatrixBlock< Scalar, Dim, Dim > >, Dune::BlockVector< Dune::FieldVector< Scalar, Dim > >, Dune::BlockVector< Dune::FieldVector< Scalar, Dim > > > OpBSeq
Definition: PreconditionerFactory_impl.hpp:203
Dune::OverlappingSchwarzOperator< Dune::BCRSMatrix< MatrixBlock< Scalar, Dim, Dim > >, Dune::BlockVector< Dune::FieldVector< Scalar, Dim > >, Dune::BlockVector< Dune::FieldVector< Scalar, Dim > >, CommPar > OpBPar
Definition: PreconditionerFactory_impl.hpp:229
Dune::MatrixAdapter< Dune::BCRSMatrix< Dune::FieldMatrix< Scalar, Dim, Dim > >, Dune::BlockVector< Dune::FieldVector< Scalar, Dim > >, Dune::BlockVector< Dune::FieldVector< Scalar, Dim > > > OpFSeq
Definition: PreconditionerFactory_impl.hpp:199
Dune::Amg::SequentialInformation CommSeq
Definition: PreconditionerFactory_impl.hpp:194
Dune::OverlappingSchwarzOperator< Dune::BCRSMatrix< Dune::FieldMatrix< Scalar, Dim, Dim > >, Dune::BlockVector< Dune::FieldVector< Scalar, Dim > >, Dune::BlockVector< Dune::FieldVector< Scalar, Dim > >, CommPar > OpFPar
Definition: PreconditionerFactory_impl.hpp:223
Dune::OwnerOverlapCopyCommunication< int, int > CommPar
Definition: PreconditionerFactory_impl.hpp:217
static void add()
Definition: StandardPreconditioners_mpi.hpp:129