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
31{
32namespace Accelerator
33{
34
36template <unsigned int block_size>
37class openclSolverBackend : public BdaSolver<block_size>
38{
40
41 using Base::N;
42 using Base::Nb;
43 using Base::nnz;
44 using Base::nnzb;
45 using Base::verbosity;
46 using Base::platformID;
47 using Base::deviceID;
48 using Base::maxit;
49 using Base::tolerance;
51
52private:
53 double *h_b = nullptr; // b vector, on host
54 std::vector<double> vals_contiguous; // only used if COPY_ROW_BY_ROW is true in openclSolverBackend.cpp
55
56 // OpenCL variables must be reusable, they are initialized in initialize()
57 cl::Buffer d_Avals, d_Acols, d_Arows; // matrix in BSR format on GPU
58 cl::Buffer d_x, d_b, d_rb, d_r, d_rw, d_p; // vectors, used during linear solve
59 cl::Buffer d_pw, d_s, d_t, d_v; // vectors, used during linear solve
60 cl::Buffer d_tmp; // used as tmp GPU buffer for dot() and norm()
61
62 std::vector<cl::Device> devices;
63
64 bool useJacMatrix = false;
65
66 std::unique_ptr<Preconditioner<block_size> > prec;
67 // can perform blocked ILU0 and AMG on pressure component
68 bool is_root; // allow for nested solvers, the root solver is called by BdaBridge
69 bool analysis_done = false;
70 std::shared_ptr<BlockedMatrix> mat = nullptr; // original matrix
71 std::shared_ptr<BlockedMatrix> jacMat = nullptr; // matrix for preconditioner
72 bool opencl_ilu_parallel; // parallelize ILU operations (with level_scheduling)
73 std::vector<cl::Event> events;
74 cl_int err;
75
79 void gpu_pbicgstab(WellContributions& wellContribs, BdaResult& res);
80
84 void initialize(std::shared_ptr<BlockedMatrix> matrix, std::shared_ptr<BlockedMatrix> jacMatrix);
85
87 void copy_system_to_gpu();
88
92 void update_system(double *vals, double *b);
93
95 void update_system_on_gpu();
96
99 bool analyze_matrix();
100
103 bool create_preconditioner();
104
109 void solve_system(WellContributions &wellContribs, BdaResult &res);
110
111public:
112 std::shared_ptr<cl::Context> context;
113 std::shared_ptr<cl::CommandQueue> queue;
114
124 openclSolverBackend(int linear_solver_verbosity, int maxit, double tolerance, unsigned int platformID, unsigned int deviceID,
125 bool opencl_ilu_parallel, std::string linsolver);
126
128 openclSolverBackend(int linear_solver_verbosity, int maxit, double tolerance, bool opencl_ilu_parallel);
129
137 SolverStatus solve_system(std::shared_ptr<BlockedMatrix> matrix, double *b,
138 std::shared_ptr<BlockedMatrix> jacMatrix, WellContributions& wellContribs, BdaResult &res) override;
139
142 // SolverStatus solve_system(BdaResult &res);
143
146 void get_result(double *x) override;
147
152 void setOpencl(std::shared_ptr<cl::Context>& context, std::shared_ptr<cl::CommandQueue>& queue);
153
154}; // end class openclSolverBackend
155
156} // namespace Accelerator
157} // namespace Opm
158
159#endif
160
161
Definition: BdaResult.hpp:31
Definition: BdaSolver.hpp:46
int nnz
Definition: BdaSolver.hpp:63
int verbosity
Definition: BdaSolver.hpp:56
int maxit
Definition: BdaSolver.hpp:58
double tolerance
Definition: BdaSolver.hpp:59
int N
Definition: BdaSolver.hpp:61
int Nb
Definition: BdaSolver.hpp:62
int nnzb
Definition: BdaSolver.hpp:64
bool initialized
Definition: BdaSolver.hpp:69
unsigned int deviceID
Definition: BdaSolver.hpp:67
unsigned int platformID
Definition: BdaSolver.hpp:66
This class implements a opencl-based ilu0-bicgstab solver on GPU.
Definition: openclSolverBackend.hpp:38
std::shared_ptr< cl::Context > context
Definition: openclSolverBackend.hpp:112
void get_result(double *x) override
openclSolverBackend(int linear_solver_verbosity, int maxit, double 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)
std::shared_ptr< cl::CommandQueue > queue
Definition: openclSolverBackend.hpp:113
SolverStatus solve_system(std::shared_ptr< BlockedMatrix > matrix, double *b, std::shared_ptr< BlockedMatrix > jacMatrix, WellContributions &wellContribs, BdaResult &res) override
openclSolverBackend(int linear_solver_verbosity, int maxit, double tolerance, bool opencl_ilu_parallel)
For the CPR coarse solver.
Definition: WellContributions.hpp:52
SolverStatus
Definition: BdaSolver.hpp:35
Definition: BlackoilPhases.hpp:27