opm-simulators
Loading...
Searching...
No Matches
ISTLSolverRuntimeOptionProxy.hpp
1/*
2 Copyright 2025 Equinor ASA
3
4 This file is part of the Open Porous Media project (OPM).
5 OPM is free software: you can redistribute it and/or modify
6 it under the terms of the GNU General Public License as published by
7 the Free Software Foundation, either version 3 of the License, or
8 (at your option) any later version.
9 OPM is distributed in the hope that it will be useful,
10 but WITHOUT ANY WARRANTY; without even the implied warranty of
11 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 GNU General Public License for more details.
13 You should have received a copy of the GNU General Public License
14 along with OPM. If not, see <http://www.gnu.org/licenses/>.
15*/
16
17#ifndef OPM_ISTLSOLVERRUNTIMEOPTIONPROXY_HEADER_INCLUDED
18#define OPM_ISTLSOLVERRUNTIMEOPTIONPROXY_HEADER_INCLUDED
19
20#include "opm/simulators/linalg/FlowLinearSolverParameters.hpp"
21#include <opm/simulators/linalg/setupPropertyTree.hpp>
22#include <opm/simulators/linalg/AbstractISTLSolver.hpp>
23#include <opm/simulators/linalg/ISTLSolver.hpp>
24#if COMPILE_GPU_BRIDGE
25#include <opm/simulators/linalg/ISTLSolverGpuBridge.hpp>
26#endif
27
28#if USE_HIP
29#include <opm/simulators/linalg/gpuistl_hip/ISTLSolverGPUISTL.hpp>
30#elif HAVE_CUDA
31#include <opm/simulators/linalg/gpuistl/ISTLSolverGPUISTL.hpp>
32#endif
33
34#include <fmt/format.h>
35
36namespace Opm
37{
38
42template <class TypeTag>
43class ISTLSolverRuntimeOptionProxy : public AbstractISTLSolver<GetPropType<TypeTag, Properties::SparseMatrixAdapter>,
44 GetPropType<TypeTag, Properties::GlobalEqVector>>
45{
46public:
50 using Matrix = typename SparseMatrixAdapter::IstlMatrix;
51
52#if HAVE_MPI
53 using CommunicationType = Dune::OwnerOverlapCopyCommunication<int, int>;
54#else
55 using CommunicationType = Dune::Communication<int>;
56#endif
57
58
59 static void registerParameters()
60 {
61 FlowLinearSolverParameters::registerParameters();
62 }
63
71 ISTLSolverRuntimeOptionProxy(const Simulator& simulator,
72 const FlowLinearSolverParameters& parameters,
73 bool forceSerial = false)
74 {
75 createSolver(simulator, parameters, forceSerial);
76 }
77
80 explicit ISTLSolverRuntimeOptionProxy(const Simulator& simulator)
81 {
82 createSolver(simulator);
83 }
84
85
86 void eraseMatrix() override
87 {
88 istlSolver_->eraseMatrix();
89 }
90
91 void setActiveSolver(int num) override
92 {
93 istlSolver_->setActiveSolver(num);
94 }
95
96 int numAvailableSolvers() const override
97 {
98 return istlSolver_->numAvailableSolvers();
99 }
100
101 void prepare(const SparseMatrixAdapter& M, Vector& b) override
102 {
103 istlSolver_->prepare(M, b);
104 }
105
106 void prepare(const Matrix& M, Vector& b) override
107 {
108 istlSolver_->prepare(M, b);
109 }
110
111 void setResidual(Vector& b) override
112 {
113 istlSolver_->setResidual(b);
114 }
115
116 void getResidual(Vector& b) const override
117 {
118 istlSolver_->getResidual(b);
119 }
120
121 void setMatrix(const SparseMatrixAdapter& M) override
122 {
123 istlSolver_->setMatrix(M);
124 }
125
126 bool solve(Vector& x) override
127 {
128 return istlSolver_->solve(x);
129 }
130
131 int iterations() const override
132 {
133 return istlSolver_->iterations();
134 }
135
136 const CommunicationType* comm() const override
137 {
138 return istlSolver_->comm();
139 }
140
141 int getSolveCount() const override
142 {
143 return istlSolver_->getSolveCount();
144 }
145
146private:
147 std::unique_ptr<AbstractISTLSolver<SparseMatrixAdapter, Vector>> istlSolver_;
148
149
150 template <class... Args>
151 void createSolver(const Simulator& simulator, Args&&... args)
152 {
153 const auto backend = Parameters::linearSolverAcceleratorTypeFromCLI();
154 if (backend == Parameters::LinearSolverAcceleratorType::CPU) {
155 // Note that for now we keep the old behavior of using the bridge solver if it is available.
156#if COMPILE_GPU_BRIDGE
157 istlSolver_ = std::make_unique<ISTLSolverGpuBridge<TypeTag>>(simulator, std::forward<Args>(args)...);
158#else
159 istlSolver_ = std::make_unique<ISTLSolver<TypeTag>>(simulator, std::forward<Args>(args)...);
160#endif
161 }
162#if HAVE_CUDA
163 else if (backend == Parameters::LinearSolverAcceleratorType::GPU) {
164 istlSolver_ = std::make_unique<gpuistl::ISTLSolverGPUISTL<TypeTag>>(simulator, std::forward<Args>(args)...);
165 }
166#endif
167 else {
168 // If we reach here, it means the backend is not supported. This could be because we have added a third backend
169 // that we need to handle. A user error would be handled in the linearSolverAcceleratorTypeFromString function called above.
170 OPM_THROW(std::invalid_argument, fmt::format(fmt::runtime("Unknown backend: {}"), Parameters::toString(backend)));
171 }
172 }
173};
174} // namespace Opm
175
176#endif
Abstract interface for ISTL solvers.
Definition AbstractISTLSolver.hpp:45
ISTLSolverRuntimeOptionProxy(const Simulator &simulator, const FlowLinearSolverParameters &parameters, bool forceSerial=false)
Construct a system solver.
Definition ISTLSolverRuntimeOptionProxy.hpp:71
int numAvailableSolvers() const override
Get the number of available solvers.
Definition ISTLSolverRuntimeOptionProxy.hpp:96
ISTLSolverRuntimeOptionProxy(const Simulator &simulator)
Construct a system solver.
Definition ISTLSolverRuntimeOptionProxy.hpp:80
void setActiveSolver(int num) override
Set the active solver by its index.
Definition ISTLSolverRuntimeOptionProxy.hpp:91
int getSolveCount() const override
Get the count of how many times the solver has been called.
Definition ISTLSolverRuntimeOptionProxy.hpp:141
const CommunicationType * comm() const override
Get the communication object used by the solver.
Definition ISTLSolverRuntimeOptionProxy.hpp:136
int iterations() const override
Get the number of iterations used in the last solve.
Definition ISTLSolverRuntimeOptionProxy.hpp:131
void eraseMatrix() override
Signals that the memory for the matrix internally in the solver could be erased.
Definition ISTLSolverRuntimeOptionProxy.hpp:86
Manages the initializing and running of time dependent problems.
Definition simulator.hh:84
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition propertysystem.hh:233
This class carries all parameters for the NewtonIterationBlackoilInterleaved class.
Definition FlowLinearSolverParameters.hpp:98