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 std::string type = prm.get<std::string>("type", "paroverilu0");
95 // We use lower case as the internal canonical representation of solver names
96 std::transform(type.begin(), type.end(), type.begin(), ::tolower);
97 auto it = creators_.find(type);
98 if (it == creators_.end()) {
99 std::ostringstream msg;
100 msg << "Preconditioner type " << type << " is not registered in the factory. Available types are: ";
101 for (const auto& prec : creators_) {
102 msg << prec.first << ' ';
103 }
104 msg << std::endl;
105 OPM_THROW(std::invalid_argument, msg.str());
106 }
107 return it->second(op, prm, weightsCalculator, pressureIndex);
108}
109
110template <class Operator, class Comm>
112PreconditionerFactory<Operator, Comm>::doCreate(const Operator& op,
113 const PropertyTree& prm,
114 const std::function<Vector()> weightsCalculator,
115 std::size_t pressureIndex,
116 const Comm& comm)
117{
118 if (!defAdded_) {
120 defAdded_ = true;
121 }
122 std::string type = prm.get<std::string>("type", "paroverilu0");
123 // We use lower case as the internal canonical representation of solver names
124 std::transform(type.begin(), type.end(), type.begin(), ::tolower);
125 auto it = parallel_creators_.find(type);
126 if (it == parallel_creators_.end()) {
127 std::ostringstream msg;
128 msg << "Parallel preconditioner type " << type << " is not registered in the factory. Available types are: ";
129 for (const auto& prec : parallel_creators_) {
130 msg << prec.first << ' ';
131 }
132 msg << std::endl;
133 OPM_THROW(std::invalid_argument, msg.str());
134 }
135 return it->second(op, prm, weightsCalculator, pressureIndex, comm);
136}
137
138template <class Operator, class Comm>
139void
140PreconditionerFactory<Operator, Comm>::doAddCreator(const std::string& type, Creator c)
141{
142 creators_[type] = c;
143}
144
145template <class Operator, class Comm>
146void
147PreconditionerFactory<Operator, Comm>::doAddCreator(const std::string& type, ParCreator c)
148{
149 parallel_creators_[type] = c;
150}
151
152template <class Operator, class Comm>
155 const PropertyTree& prm,
156 const std::function<Vector()>& weightsCalculator,
157 std::size_t pressureIndex)
158{
159 return instance().doCreate(op, prm, weightsCalculator, pressureIndex);
160}
161
162template <class Operator, class Comm>
165 const PropertyTree& prm,
166 const std::function<Vector()>& weightsCalculator,
167 const Comm& comm,
168 std::size_t pressureIndex)
169{
170 return instance().doCreate(op, prm, weightsCalculator, pressureIndex, comm);
171}
172
173
174template <class Operator, class Comm>
177 const PropertyTree& prm,
178 const Comm& comm,
179 std::size_t pressureIndex)
180{
181 return instance().doCreate(op, prm, std::function<Vector()>(), pressureIndex, comm);
182}
183
184template <class Operator, class Comm>
185void
187{
188 instance().doAddCreator(type, creator);
189}
190
191template <class Operator, class Comm>
192void
194{
195 instance().doAddCreator(type, creator);
196}
197
198using CommSeq = Dune::Amg::SequentialInformation;
199
200template<class Scalar, int Dim>
201using OpFSeq = Dune::MatrixAdapter<Dune::BCRSMatrix<Dune::FieldMatrix<Scalar, Dim, Dim>>,
202 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>,
203 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>>;
204template<class Scalar, int Dim>
205using OpBSeq = Dune::MatrixAdapter<Dune::BCRSMatrix<MatrixBlock<Scalar, Dim, Dim>>,
206 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>,
207 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>>;
208
209template<class Scalar, int Dim, bool overlap>
211 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>,
212 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>>;
213
214template<class Scalar, int Dim, bool overlap>
216 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>,
217 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>,
218 overlap>;
219
220#if HAVE_MPI
221using CommPar = Dune::OwnerOverlapCopyCommunication<int, int>;
222
223template<class Scalar, int Dim>
224using OpFPar = Dune::OverlappingSchwarzOperator<Dune::BCRSMatrix<Dune::FieldMatrix<Scalar, Dim, Dim>>,
225 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>,
226 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>,
227 CommPar>;
228
229template<class Scalar, int Dim>
230using OpBPar = Dune::OverlappingSchwarzOperator<Dune::BCRSMatrix<MatrixBlock<Scalar, Dim, Dim>>,
231 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>,
232 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>,
233 CommPar>;
234template<class Scalar, int Dim>
236 Dune::BlockVector<Dune::FieldVector<Scalar,Dim>>,
237 Dune::BlockVector<Dune::FieldVector<Scalar,Dim>>,
238 CommPar>;
239
240template<class Scalar, int Dim>
242 Dune::BlockVector<Dune::FieldVector<Scalar,Dim>>,
243 Dune::BlockVector<Dune::FieldVector<Scalar,Dim>>,
244 CommPar>;
245
246#define INSTANTIATE_PF_PAR(T,Dim) \
247 template class PreconditionerFactory<OpBSeq<T,Dim>, CommPar>; \
248 template class PreconditionerFactory<OpFPar<T,Dim>, CommPar>; \
249 template class PreconditionerFactory<OpBPar<T,Dim>, CommPar>; \
250 template class PreconditionerFactory<OpGLFPar<T,Dim>, CommPar>; \
251 template class PreconditionerFactory<OpGLBPar<T,Dim>, CommPar>; \
252 template class PreconditionerFactory<OpW<T,Dim, false>, CommPar>; \
253 template class PreconditionerFactory<OpWG<T,Dim, true>, CommPar>; \
254 template class PreconditionerFactory<OpBPar<T,Dim>, CommSeq>; \
255 template class PreconditionerFactory<OpGLBPar<T,Dim>, CommSeq>;
256#endif
257
258#define INSTANTIATE_PF_SEQ(T,Dim) \
259 template class PreconditionerFactory<OpFSeq<T,Dim>, CommSeq>; \
260 template class PreconditionerFactory<OpBSeq<T,Dim>, CommSeq>; \
261 template class PreconditionerFactory<OpW<T,Dim, false>, CommSeq>; \
262 template class PreconditionerFactory<OpWG<T,Dim, true>, CommSeq>;
263
264#if HAVE_MPI
265#define INSTANTIATE_PF(T,Dim) \
266 INSTANTIATE_PF_PAR(T,Dim) \
267 INSTANTIATE_PF_SEQ(T,Dim)
268#else
269#define INSTANTIATE_PF(T,Dim) INSTANTIATE_PF_SEQ(T,Dim)
270#endif
271
272} // 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:154
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:186
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:207
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:233
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:203
Dune::Amg::SequentialInformation CommSeq
Definition: PreconditionerFactory_impl.hpp:198
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:227
Dune::OwnerOverlapCopyCommunication< int, int > CommPar
Definition: PreconditionerFactory_impl.hpp:221
static void add()
Definition: StandardPreconditioners_mpi.hpp:129