21#ifndef OPM_FLEXIBLE_SOLVER_IMPL_HEADER_INCLUDED
22#define OPM_FLEXIBLE_SOLVER_IMPL_HEADER_INCLUDED
24#include <opm/common/ErrorMacros.hpp>
25#include <opm/common/TimingMacros.hpp>
26#include <opm/simulators/linalg/matrixblock.hh>
27#include <opm/simulators/linalg/ilufirstelement.hh>
28#include <opm/simulators/linalg/FlexibleSolver.hpp>
29#include <opm/simulators/linalg/PreconditionerFactory.hpp>
30#include <opm/simulators/linalg/PropertyTree.hpp>
31#include <opm/simulators/linalg/WellOperators.hpp>
32#include <opm/simulators/linalg/PreconditionerFactoryGPUIncludeWrapper.hpp>
33#include <opm/simulators/linalg/is_gpu_operator.hpp>
35#include <dune/common/fmatrix.hh>
36#include <dune/istl/bcrsmatrix.hh>
37#include <dune/istl/solvers.hh>
38#include <dune/istl/umfpack.hh>
39#include <dune/istl/owneroverlapcopy.hh>
40#include <dune/istl/paamg/pinfo.hh>
45#include <opm/simulators/linalg/gpuistl_hip/SolverAdapter.hpp>
47#include <opm/simulators/linalg/gpuistl/SolverAdapter.hpp>
54 template <
class Operator>
58 const std::function<VectorType()>& weightsCalculator,
59 std::size_t pressureIndex)
61 init(op, Dune::Amg::SequentialInformation(), prm, weightsCalculator,
66 template <
class Operator>
72 const std::function<VectorType()>& weightsCalculator,
73 std::size_t pressureIndex)
75 init(op, comm, prm, weightsCalculator, pressureIndex);
78 template <
class Operator>
81 apply(VectorType& x, VectorType& rhs, Dune::InverseOperatorResult& res)
84 recreateDirectSolver();
86 linsolver_->apply(x, rhs, res);
89 template <
class Operator>
92 apply(VectorType& x, VectorType& rhs,
double reduction, Dune::InverseOperatorResult& res)
95 recreateDirectSolver();
97 linsolver_->apply(x, rhs, reduction, res);
101 template <
class Operator>
106 return *preconditioner_;
109 template <
class Operator>
110 Dune::SolverCategory::Category
114 return linearoperator_for_solver_->category();
118 template <
class Operator>
119 template <
class Comm>
121 FlexibleSolver<Operator>::
122 initOpPrecSp(Operator& op,
124 const std::function<VectorType()> weightsCalculator,
126 std::size_t pressureIndex)
129 linearoperator_for_solver_ = &op;
136 scalarproduct_ = Dune::createScalarProduct<VectorType, Comm>(comm, op.category());
139 template <
class Operator>
143 const Opm::PropertyTree& prm,
144 const std::function<VectorType()> weightsCalculator,
145 const Dune::Amg::SequentialInformation&,
146 std::size_t pressureIndex)
149 linearoperator_for_solver_ = &op;
152 child ? *child : Opm::PropertyTree(),
155 scalarproduct_ = std::make_shared<Dune::SeqScalarProduct<VectorType>>();
159 template <
class Operator>
160 template <
class Comm>
163 initSolver(
const Opm::PropertyTree& prm,
const Comm& comm)
165 const bool is_iorank = comm.communicator().rank() == 0;
166 const double tol = prm.
get<
double>(
"tol", 1e-2);
167 const int maxiter = prm.
get<
int>(
"maxiter", 200);
168 const int verbosity = is_iorank ? prm.
get<
int>(
"verbosity", 0) : 0;
169 const std::string solver_type = prm.
get<std::string>(
"solver",
"bicgstab");
177 if (solver_type ==
"bicgstab") {
178 linsolver_ = std::make_shared<Dune::BiCGSTABSolver<VectorType>>(*linearoperator_for_solver_,
184 }
else if (solver_type ==
"loopsolver") {
185 linsolver_ = std::make_shared<Dune::LoopSolver<VectorType>>(*linearoperator_for_solver_,
191 }
else if (solver_type ==
"gmres") {
192 int restart = prm.
get<
int>(
"restart", 15);
193 linsolver_ = std::make_shared<Dune::RestartedGMResSolver<VectorType>>(*linearoperator_for_solver_,
202 if constexpr (!Opm::is_gpu_operator_v<Operator>) {
203 if (solver_type ==
"flexgmres") {
204 int restart = prm.
get<
int>(
"restart", 15);
205 linsolver_ = std::make_shared<Dune::RestartedFlexibleGMResSolver<VectorType>>(*linearoperator_for_solver_,
212#if HAVE_SUITESPARSE_UMFPACK
213 }
else if (solver_type ==
"umfpack") {
214 if constexpr (std::is_same_v<typename VectorType::field_type,float>) {
215 OPM_THROW(std::invalid_argument,
"UMFPack cannot be used with floats");
217 using MatrixType = std::remove_const_t<std::remove_reference_t<
decltype(linearoperator_for_solver_->getmat())>>;
218 linsolver_ = std::make_shared<Dune::UMFPack<MatrixType>>(linearoperator_for_solver_->getmat(), verbosity,
false);
219 direct_solver_ =
true;
223 }
else if (solver_type ==
"gpubicgstab") {
224 linsolver_.reset(
new Opm::gpuistl::SolverAdapter<Operator, Dune::BiCGSTABSolver, VectorType>(
225 *linearoperator_for_solver_,
237 OPM_THROW(std::invalid_argument,
238 "Properties: Solver " + solver_type +
" not known.");
248 template <
class Operator>
253#if HAVE_SUITESPARSE_UMFPACK
254 if constexpr (!Opm::is_gpu_operator_v<Operator>) {
255 if constexpr (std::is_same_v<typename VectorType::field_type, float>) {
256 OPM_THROW(std::invalid_argument,
"UMFPack cannot be used with floats");
258 using MatrixType = std::remove_const_t<std::remove_reference_t<
decltype(linearoperator_for_solver_->getmat())>>;
259 linsolver_ = std::make_shared<Dune::UMFPack<MatrixType>>(linearoperator_for_solver_->getmat(), 0,
false);
263 OPM_THROW(std::logic_error,
"Direct solver specified, but the FlexibleSolver class was not compiled with SuiteSparse support.");
270 template <
class Operator>
271 template <
class Comm>
276 const Opm::PropertyTree& prm,
277 const std::function<VectorType()> weightsCalculator,
278 std::size_t pressureIndex)
280 initOpPrecSp(op, prm, weightsCalculator, comm, pressureIndex);
281 initSolver(prm, comm);
290template<
class Scalar,
int N>
291using BV = Dune::BlockVector<Dune::FieldVector<Scalar, N>>;
292template<
class Scalar,
int N>
293using OBM = Dune::BCRSMatrix<Opm::MatrixBlock<Scalar, N, N>>;
296template<
class Scalar,
int N>
297using SeqOpM = Dune::MatrixAdapter<OBM<Scalar,N>, BV<Scalar,N>, BV<Scalar,N>>;
298template<
class Scalar,
int N>
304using Comm = Dune::OwnerOverlapCopyCommunication<int, int>;
305template<
class Scalar,
int N>
307template<
class Scalar,
int N>
309template<
class Scalar,
int N>
310using ParOpD = Dune::OverlappingSchwarzOperator<OBM<Scalar,N>, BV<Scalar,N>, BV<Scalar,N>, Comm>;
316#define INSTANTIATE_FLEXIBLESOLVER_OP(...) \
317 template class Dune::FlexibleSolver<__VA_ARGS__>; \
318 template Dune::FlexibleSolver<__VA_ARGS__>:: \
319 FlexibleSolver(__VA_ARGS__& op, \
321 const Opm::PropertyTree& prm, \
322 const std::function<typename __VA_ARGS__::domain_type()>& weightsCalculator, \
323 std::size_t pressureIndex);
325#define INSTANTIATE_FLEXIBLESOLVER(T,N) \
326 INSTANTIATE_FLEXIBLESOLVER_OP(SeqOpM<T,N>); \
327 INSTANTIATE_FLEXIBLESOLVER_OP(SeqOpW<T,N>); \
328 INSTANTIATE_FLEXIBLESOLVER_OP(ParOpM<T,N>); \
329 INSTANTIATE_FLEXIBLESOLVER_OP(ParOpW<T,N>); \
330 INSTANTIATE_FLEXIBLESOLVER_OP(ParOpD<T,N>);
334#define INSTANTIATE_FLEXIBLESOLVER_OP(...) \
335 template class Dune::FlexibleSolver<__VA_ARGS__>;
337#define INSTANTIATE_FLEXIBLESOLVER(T,N) \
338 INSTANTIATE_FLEXIBLESOLVER_OP(SeqOpM<T,N>); \
339 INSTANTIATE_FLEXIBLESOLVER_OP(SeqOpW<T,N>);
A solver class that encapsulates all needed objects for a linear solver (operator,...
Definition FlexibleSolver.hpp:45
FlexibleSolver(Operator &op, const Opm::PropertyTree &prm, const std::function< VectorType()> &weightsCalculator, std::size_t pressureIndex)
Create a sequential solver.
Definition FlexibleSolver_impl.hpp:56
Dune::PreconditionerWithUpdate< VectorType, VectorType > AbstractPrecondType
Base class type of the contained preconditioner.
Definition FlexibleSolver.hpp:50
AbstractPrecondType & preconditioner()
Access the contained preconditioner.
Definition FlexibleSolver_impl.hpp:104
Dune linear operator that assumes ghost rows are ordered after interior rows.
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())
Create a new serial preconditioner and return a pointer to it.
Definition PreconditionerFactory_impl.hpp:154
Hierarchical collection of key/value pairs.
Definition PropertyTree.hpp:39
std::optional< PropertyTree > get_child_optional(const std::string &key) const
Retrieve copy of sub tree rooted at node.
Definition PropertyTree.cpp:90
T get(const std::string &key) const
Retrieve property value given hierarchical property key.
Definition PropertyTree.cpp:59
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