Go to the documentation of this file.
23#include <opm/common/ErrorMacros.hpp>
24#include <opm/common/TimingMacros.hpp>
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>
69template < class Operator, class Comm>
70PreconditionerFactory<Operator, Comm>::PreconditionerFactory()
75template < class Operator, class Comm>
76PreconditionerFactory<Operator, Comm>&
77PreconditionerFactory<Operator, Comm>::instance()
79 static PreconditionerFactory singleton;
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)
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 << ' ';
103 OPM_THROW(std::invalid_argument, msg.str());
105 return it->second(op, prm, weightsCalculator, pressureIndex);
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,
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 << ' ';
129 OPM_THROW(std::invalid_argument, msg.str());
131 return it->second(op, prm, weightsCalculator, pressureIndex, comm);
134template < class Operator, class Comm>
136PreconditionerFactory<Operator, Comm>::doAddCreator( const std::string& type, Creator c)
141template < class Operator, class Comm>
143PreconditionerFactory<Operator, Comm>::doAddCreator( const std::string& type, ParCreator c)
145 parallel_creators_[type] = c;
148template < class Operator, class Comm>
152 const std::function< Vector()>& weightsCalculator,
153 std::size_t pressureIndex)
155 return instance().doCreate(op, prm, weightsCalculator, pressureIndex);
158template < class Operator, class Comm>
162 const std::function< Vector()>& weightsCalculator,
164 std::size_t pressureIndex)
166 return instance().doCreate(op, prm, weightsCalculator, pressureIndex, comm);
170template < class Operator, class Comm>
175 std::size_t pressureIndex)
177 return instance().doCreate(op, prm, std::function< Vector()>(), pressureIndex, comm);
180template < class Operator, class Comm>
184 instance().doAddCreator(type, creator);
187template < class Operator, class Comm>
191 instance().doAddCreator(type, creator);
194using CommSeq = Dune::Amg::SequentialInformation;
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>>>;
205template< class Scalar, int Dim, bool overlap>
207 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>,
208 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>>;
210template< class Scalar, int Dim, bool overlap>
212 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>,
213 Dune::BlockVector<Dune::FieldVector<Scalar, Dim>>,
217using CommPar = Dune::OwnerOverlapCopyCommunication<int, int>;
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>>,
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>>,
230template< class Scalar, int Dim>
232 Dune::BlockVector<Dune::FieldVector<Scalar,Dim>>,
233 Dune::BlockVector<Dune::FieldVector<Scalar,Dim>>,
236template< class Scalar, int Dim>
238 Dune::BlockVector<Dune::FieldVector<Scalar,Dim>>,
239 Dune::BlockVector<Dune::FieldVector<Scalar,Dim>>,
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>;
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>;
261#define INSTANTIATE_PF(T,Dim) \
262 INSTANTIATE_PF_PAR(T,Dim) \
263 INSTANTIATE_PF_SEQ(T,Dim)
265#define INSTANTIATE_PF(T,Dim) INSTANTIATE_PF_SEQ(T,Dim)
Dune::OwnerOverlapCopyCommunication< int, int > Comm Definition: FlexibleSolver_impl.hpp:304
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
|