SolverAdapter.hpp
Go to the documentation of this file.
1/*
2 Copyright 2022-2023 SINTEF AS
3
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#ifndef OPM_SOLVERADAPTER_HPP
20#define OPM_SOLVERADAPTER_HPP
21
22#include <memory>
23
24#include <dune/istl/bcrsmatrix.hh>
25#include <dune/istl/operators.hh>
26#include <dune/istl/paamg/pinfo.hh>
27#include <dune/istl/schwarz.hh>
28#include <dune/istl/solver.hh>
29#include <opm/common/ErrorMacros.hpp>
36
37#ifdef OPEN_MPI
38#if OPEN_MPI
39#include "mpi-ext.h"
40#endif
41#endif
42
43namespace Opm::gpuistl
44{
50template <class Operator, template <class> class UnderlyingSolver, class X>
51class SolverAdapter : public Dune::IterativeSolver<X, X>
52{
53public:
54 using typename Dune::IterativeSolver<X, X>::domain_type;
55 using typename Dune::IterativeSolver<X, X>::range_type;
56 using typename Dune::IterativeSolver<X, X>::field_type;
57 using typename Dune::IterativeSolver<X, X>::real_type;
58 using typename Dune::IterativeSolver<X, X>::scalar_real_type;
59 static constexpr auto block_size = domain_type::block_type::dimension;
60 using XGPU = Opm::gpuistl::GpuVector<real_type>;
61
62
73 template<class Comm>
74 SolverAdapter(Operator& op,
75 Dune::ScalarProduct<X>& sp,
76 std::shared_ptr<Dune::Preconditioner<X, X>> prec,
77 scalar_real_type reduction,
78 int maxit,
79 int verbose,
80 const Comm& comm)
81 : Dune::IterativeSolver<X, X>(op, sp, *prec, reduction, maxit, verbose)
82 , m_opOnCPUWithMatrix(op)
83 , m_matrix(GpuSparseMatrix<real_type>::fromMatrix(op.getmat()))
84 , m_underlyingSolver(constructSolver(prec, reduction, maxit, verbose, comm))
85 {
86 OPM_ERROR_IF(
88 "Currently we only support operators of type MatrixAdapter in the CUDA/HIP solver. "
89 "Use --matrix-add-well-contributions=true. "
90 "Using WellModelMatrixAdapter with SolverAdapter is not well-defined so it will not work well, or at all.");
91 }
92
93 virtual void apply(X& x, X& b, double reduction, Dune::InverseOperatorResult& res) override
94 {
95 // TODO: Can we do this without reimplementing the other function?
96 // TODO: [perf] Do we need to update the matrix every time? Probably yes
97 m_matrix.updateNonzeroValues(m_opOnCPUWithMatrix.getmat());
98
99 if (!m_inputBuffer) {
100 m_inputBuffer.reset(new XGPU(b.dim()));
101 m_outputBuffer.reset(new XGPU(x.dim()));
102 }
103
104 m_inputBuffer->copyFromHost(b);
105 // TODO: [perf] do we need to copy x here?
106 m_outputBuffer->copyFromHost(x);
107
108 m_underlyingSolver.apply(*m_outputBuffer, *m_inputBuffer, reduction, res);
109
110 // TODO: [perf] do we need to copy b here?
111 m_inputBuffer->copyToHost(b);
112 m_outputBuffer->copyToHost(x);
113 }
114 virtual void apply(X& x, X& b, Dune::InverseOperatorResult& res) override
115 {
116 // TODO: [perf] Do we need to update the matrix every time? Probably yes
117 m_matrix.updateNonzeroValues(m_opOnCPUWithMatrix.getmat());
118
119 if (!m_inputBuffer) {
120 m_inputBuffer.reset(new XGPU(b.dim()));
121 m_outputBuffer.reset(new XGPU(x.dim()));
122 }
123
124 m_inputBuffer->copyFromHost(b);
125 // TODO: [perf] do we need to copy x here?
126 m_outputBuffer->copyFromHost(x);
127
128 m_underlyingSolver.apply(*m_outputBuffer, *m_inputBuffer, res);
129
130 // TODO: [perf] do we need to copy b here?
131 m_inputBuffer->copyToHost(b);
132 m_outputBuffer->copyToHost(x);
133 }
134
135private:
136 Operator& m_opOnCPUWithMatrix;
138
139 UnderlyingSolver<XGPU> m_underlyingSolver;
140
141
142 // TODO: Use a std::forward
143 // This is the MPI parallel case (general communication object)
144 template <class Comm>
145 UnderlyingSolver<XGPU> constructSolver(std::shared_ptr<Dune::Preconditioner<X, X>> prec,
146 scalar_real_type reduction,
147 int maxit,
148 int verbose,
149 const Comm& communication)
150 {
151 // TODO: See the below TODO over the definition of precHolder in the other overload of constructSolver
152 // TODO: We are currently double wrapping preconditioners in the preconditioner factory to be extra
153 // compatible with CPU. Probably a cleaner way eventually would be to do more modifications to the
154 // flexible solver to accomodate the pure GPU better.
155 auto precAsHolder = std::dynamic_pointer_cast<PreconditionerHolder<X, X>>(prec);
156 if (!precAsHolder) {
157 OPM_THROW(std::invalid_argument,
158 "The preconditioner needs to be a CUDA preconditioner (eg. GPUDILU) wrapped in a "
159 "Opm::gpuistl::PreconditionerAdapter wrapped in a "
160 "Opm::gpuistl::GpuBlockPreconditioner. If you are unsure what this means, set "
161 "preconditioner to 'GPUDILU'");
162 }
163
164 auto preconditionerAdapter = precAsHolder->getUnderlyingPreconditioner();
165 auto preconditionerAdapterAsHolder
166 = std::dynamic_pointer_cast<PreconditionerHolder<XGPU, XGPU>>(preconditionerAdapter);
167 if (!preconditionerAdapterAsHolder) {
168 OPM_THROW(std::invalid_argument,
169 "The preconditioner needs to be a CUDA preconditioner (eg. GPUDILU) wrapped in a "
170 "Opm::gpuistl::PreconditionerAdapter wrapped in a "
171 "Opm::gpuistl::GpuBlockPreconditioner. If you are unsure what this means, set "
172 "preconditioner to 'GPUDILU'");
173 }
174 // We need to get the underlying preconditioner:
175 auto preconditionerReallyOnGPU = preconditionerAdapterAsHolder->getUnderlyingPreconditioner();
176
177 // Temporary solution use the GPU Direct communication solely based on these prepcrosessor statements
178 bool mpiSupportsCudaAwareAtCompileTime = false;
179 bool mpiSupportsCudaAwareAtRunTime = false;
180
181#if defined(MPIX_CUDA_AWARE_SUPPORT) && MPIX_CUDA_AWARE_SUPPORT
182 mpiSupportsCudaAwareAtCompileTime = true;
183#endif /* MPIX_CUDA_AWARE_SUPPORT */
184
185#if defined(MPIX_CUDA_AWARE_SUPPORT)
186 if (1 == MPIX_Query_cuda_support()) {
187 mpiSupportsCudaAwareAtRunTime = true;
188 }
189#endif /* MPIX_CUDA_AWARE_SUPPORT */
190
191
192 // TODO add typename Operator communication type as a named type with using
193 std::shared_ptr<Opm::gpuistl::GPUSender<real_type, Comm>> gpuComm;
194 if (mpiSupportsCudaAwareAtCompileTime && mpiSupportsCudaAwareAtRunTime) {
195 gpuComm = std::make_shared<
197 communication);
198 } else {
199 gpuComm = std::make_shared<
201 communication);
202 }
203
204 using CudaCommunication = GpuOwnerOverlapCopy<real_type, block_size, Comm>;
205 using SchwarzOperator
206 = Dune::OverlappingSchwarzOperator<GpuSparseMatrix<real_type>, XGPU, XGPU, CudaCommunication>;
207 auto cudaCommunication = std::make_shared<CudaCommunication>(gpuComm);
208
209 auto mpiPreconditioner = std::make_shared<GpuBlockPreconditioner<XGPU, XGPU, CudaCommunication>>(
210 preconditionerReallyOnGPU, cudaCommunication);
211
212 auto scalarProduct = std::make_shared<Dune::ParallelScalarProduct<XGPU, CudaCommunication>>(
213 cudaCommunication, m_opOnCPUWithMatrix.category());
214
215
216 // NOTE: Ownership of cudaCommunication is handled by mpiPreconditioner. However, just to make sure we
217 // remember this, we add this check. So remember that we hold one count in this scope, and one in the
218 // GpuBlockPreconditioner instance. We accommodate for the fact that it could be passed around in
219 // GpuBlockPreconditioner, hence we do not test for != 2
220 OPM_ERROR_IF(cudaCommunication.use_count() < 2, "Internal error. Shared pointer not owned properly.");
221 auto overlappingCudaOperator = std::make_shared<SchwarzOperator>(m_matrix, *cudaCommunication);
222
223 return UnderlyingSolver<XGPU>(
224 overlappingCudaOperator, scalarProduct, mpiPreconditioner, reduction, maxit, verbose);
225 }
226
227 // This is the serial case (specific overload for Dune::Amg::SequentialInformation)
228 UnderlyingSolver<XGPU> constructSolver(std::shared_ptr<Dune::Preconditioner<X, X>> prec,
229 scalar_real_type reduction,
230 int maxit,
231 int verbose,
232 [[maybe_unused]] const Dune::Amg::SequentialInformation& communication)
233 {
234 // Dune::Amg::SequentialInformation is the serial case
235 return constructSolver(prec, reduction, maxit, verbose);
236 }
237
238 // TODO: Use a std::forward
239 UnderlyingSolver<XGPU> constructSolver(std::shared_ptr<Dune::Preconditioner<X, X>> prec,
240 scalar_real_type reduction,
241 int maxit,
242 int verbose)
243 {
244 // TODO: Fix the reliance on casting here. This is a code smell to a certain degree, and assumes
245 // a certain setup beforehand. The only reason we do it this way is to minimize edits to the
246 // flexible solver. We could design it differently, but keep this for the time being until
247 // we figure out how we want to GPU-ify the rest of the system.
248 auto precAsHolder = std::dynamic_pointer_cast<PreconditionerHolder<XGPU, XGPU>>(prec);
249 if (!precAsHolder) {
250 OPM_THROW(std::invalid_argument,
251 "The preconditioner needs to be a CUDA preconditioner wrapped in a "
252 "Opm::gpuistl::PreconditionerHolder (eg. GPUDILU).");
253 }
254 auto preconditionerOnGPU = precAsHolder->getUnderlyingPreconditioner();
255
256 auto matrixOperator = std::make_shared<Dune::MatrixAdapter<GpuSparseMatrix<real_type>, XGPU, XGPU>>(m_matrix);
257 auto scalarProduct = std::make_shared<Dune::SeqScalarProduct<XGPU>>();
258 return UnderlyingSolver<XGPU>(matrixOperator, scalarProduct, preconditionerOnGPU, reduction, maxit, verbose);
259 }
260
261 std::unique_ptr<XGPU> m_inputBuffer;
262 std::unique_ptr<XGPU> m_outputBuffer;
263};
264} // namespace Opm::gpuistl
265
266#endif
Dune::OwnerOverlapCopyCommunication< int, int > Comm
Definition: FlexibleSolver_impl.hpp:285
Derived class of GPUSender that handles MPI made with CUDA aware MPI The copOwnerToAll function uses ...
Definition: GpuOwnerOverlapCopy.hpp:169
Derived class of GPUSender that handles MPI calls that should NOT use GPU direct communicatoin The im...
Definition: GpuOwnerOverlapCopy.hpp:114
The GpuSparseMatrix class simple wrapper class for a CuSparse matrix.
Definition: GpuSparseMatrix.hpp:48
void updateNonzeroValues(const MatrixType &matrix, bool copyNonZeroElementsDirectly=false)
updateNonzeroValues updates the non-zero values by using the non-zero values of the supplied matrix
Wraps a CUDA solver to work with CPU data.
Definition: SolverAdapter.hpp:52
virtual void apply(X &x, X &b, Dune::InverseOperatorResult &res) override
Definition: SolverAdapter.hpp:114
virtual void apply(X &x, X &b, double reduction, Dune::InverseOperatorResult &res) override
Definition: SolverAdapter.hpp:93
Opm::gpuistl::GpuVector< real_type > XGPU
Definition: SolverAdapter.hpp:60
SolverAdapter(Operator &op, Dune::ScalarProduct< X > &sp, std::shared_ptr< Dune::Preconditioner< X, X > > prec, scalar_real_type reduction, int maxit, int verbose, const Comm &comm)
constructor
Definition: SolverAdapter.hpp:74
static constexpr auto block_size
Definition: SolverAdapter.hpp:59
The is_a_well_operator class tries to guess if the operator is a well type operator.
Definition: has_function.hpp:114
Definition: fvbaseprimaryvariables.hh:141
Definition: autotuner.hpp:29
Dune::InverseOperatorResult InverseOperatorResult
Definition: BdaBridge.hpp:32