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#if HAVE_AVX2_EXTENSION
36#include <opm/simulators/linalg/mixed/wrapper.hpp>
39#include <dune/common/fmatrix.hh>
40#include <dune/istl/bcrsmatrix.hh>
41#include <dune/istl/solvers.hh>
42#include <dune/istl/umfpack.hh>
43#include <dune/istl/owneroverlapcopy.hh>
44#include <dune/istl/paamg/pinfo.hh>
49#include <opm/simulators/linalg/gpuistl_hip/SolverAdapter.hpp>
51#include <opm/simulators/linalg/gpuistl/SolverAdapter.hpp>
58 template <
class Operator>
62 const std::function<VectorType()>& weightsCalculator,
63 std::size_t pressureIndex)
65 init(op, Dune::Amg::SequentialInformation(), prm, weightsCalculator,
70 template <
class Operator>
76 const std::function<VectorType()>& weightsCalculator,
77 std::size_t pressureIndex)
79 init(op, comm, prm, weightsCalculator, pressureIndex);
82 template <
class Operator>
85 apply(VectorType& x, VectorType& rhs, Dune::InverseOperatorResult& res)
88 recreateDirectSolver();
90 linsolver_->apply(x, rhs, res);
93 template <
class Operator>
96 apply(VectorType& x, VectorType& rhs,
double reduction, Dune::InverseOperatorResult& res)
99 recreateDirectSolver();
101 linsolver_->apply(x, rhs, reduction, res);
105 template <
class Operator>
110 return *preconditioner_;
113 template <
class Operator>
114 Dune::SolverCategory::Category
118 return linearoperator_for_solver_->category();
122 template <
class Operator>
123 template <
class Comm>
125 FlexibleSolver<Operator>::
126 initOpPrecSp(Operator& op,
128 const std::function<VectorType()> weightsCalculator,
130 std::size_t pressureIndex)
133 linearoperator_for_solver_ = &op;
140 scalarproduct_ = Dune::createScalarProduct<VectorType, Comm>(comm, op.category());
143 template <
class Operator>
147 const Opm::PropertyTree& prm,
148 const std::function<VectorType()> weightsCalculator,
149 const Dune::Amg::SequentialInformation&,
150 std::size_t pressureIndex)
153 linearoperator_for_solver_ = &op;
156 child ? *child : Opm::PropertyTree(),
159 scalarproduct_ = std::make_shared<Dune::SeqScalarProduct<VectorType>>();
163 template <
class Operator>
164 template <
class Comm>
167 initSolver(
const Opm::PropertyTree& prm,
const Comm& comm)
169 const bool is_iorank = comm.communicator().rank() == 0;
170 const double tol = prm.
get<
double>(
"tol", 1e-2);
171 const int maxiter = prm.
get<
int>(
"maxiter", 200);
172 const int verbosity = is_iorank ? prm.
get<
int>(
"verbosity", 0) : 0;
173 const std::string solver_type = prm.
get<std::string>(
"solver",
"bicgstab");
181 if (solver_type ==
"bicgstab") {
182 linsolver_ = std::make_shared<Dune::BiCGSTABSolver<VectorType>>(*linearoperator_for_solver_,
188#if HAVE_AVX2_EXTENSION
189 }
else if (solver_type ==
"mixed-bicgstab") {
190 if constexpr (Opm::is_gpu_operator_v<Operator>) {
191 OPM_THROW(std::invalid_argument,
"mixed-bicgstab solver not supported for GPU operatorsg");
192 }
else if constexpr (std::is_same_v<typename VectorType::field_type, float>){
193 OPM_THROW(std::invalid_argument,
"mixed-bicgstab solver not supported for single precision.");
195 const std::string prec_type = prm.
get<std::string>(
"preconditioner.type",
"error");
196 bool use_mixed_dilu= (prec_type==
"mixed-dilu");
197 using MatrixType =
decltype(linearoperator_for_solver_->getmat());
198 linsolver_ = std::make_shared<Dune::MixedSolver<VectorType,MatrixType>>(
199 linearoperator_for_solver_->getmat(),
206 }
else if (solver_type ==
"loopsolver") {
207 linsolver_ = std::make_shared<Dune::LoopSolver<VectorType>>(*linearoperator_for_solver_,
213 }
else if (solver_type ==
"gmres") {
214 int restart = prm.
get<
int>(
"restart", 15);
215 linsolver_ = std::make_shared<Dune::RestartedGMResSolver<VectorType>>(*linearoperator_for_solver_,
223 if constexpr (!Opm::is_gpu_operator_v<Operator>) {
224 if (solver_type ==
"flexgmres") {
225 int restart = prm.
get<
int>(
"restart", 15);
226 linsolver_ = std::make_shared<Dune::RestartedFlexibleGMResSolver<VectorType>>(*linearoperator_for_solver_,
233#if HAVE_SUITESPARSE_UMFPACK
234 }
else if (solver_type ==
"umfpack") {
235 if constexpr (std::is_same_v<typename VectorType::field_type,float>) {
236 OPM_THROW(std::invalid_argument,
"UMFPack cannot be used with floats");
238 using MatrixType = std::remove_const_t<std::remove_reference_t<
decltype(linearoperator_for_solver_->getmat())>>;
239 linsolver_ = std::make_shared<Dune::UMFPack<MatrixType>>(linearoperator_for_solver_->getmat(), verbosity,
false);
240 direct_solver_ =
true;
244 }
else if (solver_type ==
"gpubicgstab") {
245 linsolver_.reset(
new Opm::gpuistl::SolverAdapter<Operator, Dune::BiCGSTABSolver, VectorType>(
246 *linearoperator_for_solver_,
258 OPM_THROW(std::invalid_argument,
259 "Properties: Solver " + solver_type +
" not known.");
269 template <
class Operator>
274#if HAVE_SUITESPARSE_UMFPACK
275 if constexpr (!Opm::is_gpu_operator_v<Operator>) {
276 if constexpr (std::is_same_v<typename VectorType::field_type, float>) {
277 OPM_THROW(std::invalid_argument,
"UMFPack cannot be used with floats");
279 using MatrixType = std::remove_const_t<std::remove_reference_t<
decltype(linearoperator_for_solver_->getmat())>>;
280 linsolver_ = std::make_shared<Dune::UMFPack<MatrixType>>(linearoperator_for_solver_->getmat(), 0,
false);
284 OPM_THROW(std::logic_error,
"Direct solver specified, but the FlexibleSolver class was not compiled with SuiteSparse support.");
291 template <
class Operator>
292 template <
class Comm>
297 const Opm::PropertyTree& prm,
298 const std::function<VectorType()> weightsCalculator,
299 std::size_t pressureIndex)
301 initOpPrecSp(op, prm, weightsCalculator, comm, pressureIndex);
302 initSolver(prm, comm);
311template<
class Scalar,
int N>
312using BV = Dune::BlockVector<Dune::FieldVector<Scalar, N>>;
313template<
class Scalar,
int N>
314using OBM = Dune::BCRSMatrix<Opm::MatrixBlock<Scalar, N, N>>;
317template<
class Scalar,
int N>
318using SeqOpM = Dune::MatrixAdapter<OBM<Scalar,N>, BV<Scalar,N>, BV<Scalar,N>>;
319template<
class Scalar,
int N>
325using Comm = Dune::OwnerOverlapCopyCommunication<int, int>;
326template<
class Scalar,
int N>
328template<
class Scalar,
int N>
330template<
class Scalar,
int N>
331using ParOpD = Dune::OverlappingSchwarzOperator<OBM<Scalar,N>, BV<Scalar,N>, BV<Scalar,N>, Comm>;
337#define INSTANTIATE_FLEXIBLESOLVER_OP(...) \
338 template class Dune::FlexibleSolver<__VA_ARGS__>; \
339 template Dune::FlexibleSolver<__VA_ARGS__>:: \
340 FlexibleSolver(__VA_ARGS__& op, \
342 const Opm::PropertyTree& prm, \
343 const std::function<typename __VA_ARGS__::domain_type()>& weightsCalculator, \
344 std::size_t pressureIndex);
346#define INSTANTIATE_FLEXIBLESOLVER(T,N) \
347 INSTANTIATE_FLEXIBLESOLVER_OP(SeqOpM<T,N>); \
348 INSTANTIATE_FLEXIBLESOLVER_OP(SeqOpW<T,N>); \
349 INSTANTIATE_FLEXIBLESOLVER_OP(ParOpM<T,N>); \
350 INSTANTIATE_FLEXIBLESOLVER_OP(ParOpW<T,N>); \
351 INSTANTIATE_FLEXIBLESOLVER_OP(ParOpD<T,N>);
355#define INSTANTIATE_FLEXIBLESOLVER_OP(...) \
356 template class Dune::FlexibleSolver<__VA_ARGS__>;
358#define INSTANTIATE_FLEXIBLESOLVER(T,N) \
359 INSTANTIATE_FLEXIBLESOLVER_OP(SeqOpM<T,N>); \
360 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:60
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:108
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