openclSolverBackend.hpp
Go to the documentation of this file.
1/*
2 Copyright 2020 Equinor ASA
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
20#ifndef OPM_OPENCLSOLVER_BACKEND_HEADER_INCLUDED
21#define OPM_OPENCLSOLVER_BACKEND_HEADER_INCLUDED
22
27
29
30namespace Opm::Accelerator {
31
33template<class Scalar, unsigned int block_size>
34class openclSolverBackend : public GpuSolver<Scalar,block_size>
35{
37
38 using Base::N;
39 using Base::Nb;
40 using Base::nnz;
41 using Base::nnzb;
42 using Base::verbosity;
43 using Base::platformID;
44 using Base::deviceID;
45 using Base::maxit;
46 using Base::tolerance;
48
49private:
50 Scalar* h_b = nullptr; // b vector, on host
51 std::vector<Scalar> vals_contiguous; // only used if COPY_ROW_BY_ROW is true in openclSolverBackend.cpp
52
53 // OpenCL variables must be reusable, they are initialized in initialize()
54 cl::Buffer d_Avals, d_Acols, d_Arows; // matrix in BSR format on GPU
55 cl::Buffer d_x, d_b, d_rb, d_r, d_rw, d_p; // vectors, used during linear solve
56 cl::Buffer d_pw, d_s, d_t, d_v; // vectors, used during linear solve
57 cl::Buffer d_tmp; // used as tmp GPU buffer for dot() and norm()
58
59 std::vector<cl::Device> devices;
60
61 bool useJacMatrix = false;
62
63 std::unique_ptr<openclPreconditioner<Scalar,block_size>> prec;
64 // can perform blocked ILU0 and AMG on pressure component
65 bool is_root; // allow for nested solvers, the root solver is called by GpuBridge
66 bool analysis_done = false;
67 std::shared_ptr<BlockedMatrix<Scalar>> mat{}; // original matrix
68 std::shared_ptr<BlockedMatrix<Scalar>> jacMat{}; // matrix for preconditioner
69 bool opencl_ilu_parallel; // parallelize ILU operations (with level_scheduling)
70 std::vector<cl::Event> events;
71 cl_int err;
72
76 void gpu_pbicgstab(WellContributions<Scalar>& wellContribs, GpuResult& res);
77
81 void initialize(std::shared_ptr<BlockedMatrix<Scalar>> matrix,
82 std::shared_ptr<BlockedMatrix<Scalar>> jacMatrix);
83
85 void copy_system_to_gpu();
86
90 void update_system(Scalar* vals, Scalar* b);
91
93 void update_system_on_gpu();
94
97 bool analyze_matrix();
98
101 bool create_preconditioner();
102
107 void solve_system(WellContributions<Scalar>& wellContribs, GpuResult& res);
108
109public:
110 std::shared_ptr<cl::Context> context{};
111 std::shared_ptr<cl::CommandQueue> queue{};
112
122 openclSolverBackend(int linear_solver_verbosity, int maxit, Scalar tolerance,
123 unsigned int platformID, unsigned int deviceID,
124 bool opencl_ilu_parallel, std::string linsolver);
125
127 openclSolverBackend(int linear_solver_verbosity, int maxit,
128 Scalar tolerance, bool opencl_ilu_parallel);
129
138 Scalar* b,
139 std::shared_ptr<BlockedMatrix<Scalar>> jacMatrix,
140 WellContributions<Scalar>& wellContribs,
141 GpuResult& res) override;
142
145 // SolverStatus solve_system(GpuResult &res);
146
149 void get_result(Scalar* x) override;
150
155 void setOpencl(std::shared_ptr<cl::Context>& context,
156 std::shared_ptr<cl::CommandQueue>& queue);
157}; // end class openclSolverBackend
158
159} // namespace Opm::Accelerator
160
161#endif
Definition: BlockedMatrix.hpp:29
Definition: GpuResult.hpp:31
Definition: GpuSolver.hpp:46
int nnzb
Definition: GpuSolver.hpp:61
unsigned int platformID
Definition: GpuSolver.hpp:63
unsigned int deviceID
Definition: GpuSolver.hpp:64
Scalar tolerance
Definition: GpuSolver.hpp:56
bool initialized
Definition: GpuSolver.hpp:66
int maxit
Definition: GpuSolver.hpp:55
int N
Definition: GpuSolver.hpp:58
int verbosity
Definition: GpuSolver.hpp:53
int nnz
Definition: GpuSolver.hpp:60
int Nb
Definition: GpuSolver.hpp:59
This class implements a opencl-based ilu0-bicgstab solver on GPU.
Definition: openclSolverBackend.hpp:35
openclSolverBackend(int linear_solver_verbosity, int maxit, Scalar tolerance, unsigned int platformID, unsigned int deviceID, bool opencl_ilu_parallel, std::string linsolver)
void setOpencl(std::shared_ptr< cl::Context > &context, std::shared_ptr< cl::CommandQueue > &queue)
SolverStatus solve_system(std::shared_ptr< BlockedMatrix< Scalar > > matrix, Scalar *b, std::shared_ptr< BlockedMatrix< Scalar > > jacMatrix, WellContributions< Scalar > &wellContribs, GpuResult &res) override
std::shared_ptr< cl::Context > context
Definition: openclSolverBackend.hpp:110
openclSolverBackend(int linear_solver_verbosity, int maxit, Scalar tolerance, bool opencl_ilu_parallel)
For the CPR coarse solver.
void get_result(Scalar *x) override
std::shared_ptr< cl::CommandQueue > queue
Definition: openclSolverBackend.hpp:111
Definition: WellContributions.hpp:51
Definition: amgclSolverBackend.hpp:44
SolverStatus
Definition: GpuSolver.hpp:35