My Project
openclSolverBackend.hpp
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 
23 #include <opm/simulators/linalg/bda/opencl/opencl.hpp>
24 #include <opm/simulators/linalg/bda/BdaResult.hpp>
25 #include <opm/simulators/linalg/bda/BdaSolver.hpp>
26 #include <opm/simulators/linalg/bda/ILUReorder.hpp>
27 #include <opm/simulators/linalg/bda/WellContributions.hpp>
28 
29 #include <opm/simulators/linalg/bda/opencl/Preconditioner.hpp>
30 
31 namespace Opm
32 {
33 namespace Accelerator
34 {
35 
37 template <unsigned int block_size>
38 class openclSolverBackend : public BdaSolver<block_size>
39 {
41 
42  using Base::N;
43  using Base::Nb;
44  using Base::nnz;
45  using Base::nnzb;
46  using Base::verbosity;
47  using Base::platformID;
48  using Base::deviceID;
49  using Base::maxit;
50  using Base::tolerance;
51  using Base::initialized;
52 
53 private:
54  double *rb = nullptr; // reordered b vector, if the matrix is reordered, rb is newly allocated, otherwise it just points to b
55  std::vector<double> vals_contiguous; // only used if COPY_ROW_BY_ROW is true in openclSolverBackend.cpp
56 
57  // OpenCL variables must be reusable, they are initialized in initialize()
58  cl::Buffer d_Avals, d_Acols, d_Arows; // (reordered) matrix in BSR format on GPU
59  cl::Buffer d_x, d_b, d_rb, d_r, d_rw, d_p; // vectors, used during linear solve
60  cl::Buffer d_pw, d_s, d_t, d_v; // vectors, used during linear solve
61  cl::Buffer d_tmp; // used as tmp GPU buffer for dot() and norm()
62  cl::Buffer d_toOrder; // only used when reordering is used
63 
64  std::vector<cl::Device> devices;
65 
66  bool useJacMatrix = false;
67 
68  std::unique_ptr<Preconditioner<block_size> > prec;
69  // can perform blocked ILU0 and AMG on pressure component
70  bool is_root; // allow for nested solvers, the root solver is called by BdaBridge
71  int *toOrder = nullptr, *fromOrder = nullptr; // BILU0 reorders rows of the matrix via these mappings
72  bool analysis_done = false;
73  std::shared_ptr<BlockedMatrix> mat = nullptr; // original matrix
74  std::shared_ptr<BlockedMatrix> jacMat = nullptr; // matrix for preconditioner
75  BlockedMatrix *rmat = nullptr; // reordered matrix (or original if no reordering), used for spmv
76  ILUReorder opencl_ilu_reorder; // reordering strategy
77  std::vector<cl::Event> events;
78  cl_int err;
79 
84  unsigned int ceilDivision(const unsigned int A, const unsigned int B);
85 
91  double dot_w(cl::Buffer in1, cl::Buffer in2, cl::Buffer out);
92 
98  double norm_w(cl::Buffer in, cl::Buffer out);
99 
104  void axpy_w(cl::Buffer in, const double a, cl::Buffer out);
105 
109  void scale_w(cl::Buffer vec, const double a);
110 
118  void custom_w(cl::Buffer p, cl::Buffer v, cl::Buffer r, const double omega, const double beta);
119 
128  void spmv_blocked_w(cl::Buffer vals, cl::Buffer cols, cl::Buffer rows, cl::Buffer x, cl::Buffer b);
129 
133  void gpu_pbicgstab(WellContributions& wellContribs, BdaResult& res);
134 
138  void initialize(std::shared_ptr<BlockedMatrix> matrix, std::shared_ptr<BlockedMatrix> jacMatrix);
139 
141  void finalize();
142 
144  void copy_system_to_gpu();
145 
150  void update_system(double *vals, double *b, WellContributions &wellContribs);
151 
153  void update_system_on_gpu();
154 
157  bool analyze_matrix();
158 
161  bool create_preconditioner();
162 
166  void solve_system(WellContributions &wellContribs, BdaResult &res);
167 
168 public:
169  std::shared_ptr<cl::Context> context;
170  std::shared_ptr<cl::CommandQueue> queue;
171 
181  openclSolverBackend(int linear_solver_verbosity, int maxit, double tolerance, unsigned int platformID, unsigned int deviceID,
182  ILUReorder opencl_ilu_reorder, std::string linsolver);
183 
185  openclSolverBackend(int linear_solver_verbosity, int maxit, double tolerance, ILUReorder opencl_ilu_reorder);
186 
189 
197  SolverStatus solve_system(std::shared_ptr<BlockedMatrix> matrix, double *b,
198  std::shared_ptr<BlockedMatrix> jacMatrix, WellContributions& wellContribs, BdaResult &res) override;
199 
202  // SolverStatus solve_system(BdaResult &res);
203 
206  void get_result(double *x) override;
207 
212  void setOpencl(std::shared_ptr<cl::Context>& context, std::shared_ptr<cl::CommandQueue>& queue);
213 
214 }; // end class openclSolverBackend
215 
216 } // namespace Accelerator
217 } // namespace Opm
218 
219 #endif
220 
221 
This class is based on InverseOperatorResult struct from dune/istl/solver.hh It is needed to prevent ...
Definition: BdaResult.hpp:31
This class serves to simplify choosing between different backend solvers, such as cusparseSolver and ...
Definition: BdaSolver.hpp:45
This struct resembles a blocked csr matrix, like Dune::BCRSMatrix.
Definition: BlockedMatrix.hpp:37
This class implements a opencl-based ilu0-bicgstab solver on GPU.
Definition: openclSolverBackend.hpp:39
~openclSolverBackend()
Destroy a openclSolver, and free memory.
Definition: openclSolverBackend.cpp:236
void get_result(double *x) override
Solve scalar linear system, for example a coarse system of an AMG preconditioner Data is already on t...
Definition: openclSolverBackend.cpp:670
openclSolverBackend(int linear_solver_verbosity, int maxit, double tolerance, unsigned int platformID, unsigned int deviceID, ILUReorder opencl_ilu_reorder, std::string linsolver)
Construct a openclSolver.
Definition: openclSolverBackend.cpp:50
void setOpencl(std::shared_ptr< cl::Context > &context, std::shared_ptr< cl::CommandQueue > &queue)
Set OpenCL objects This class either creates them based on platformID and deviceID or receives them t...
Definition: openclSolverBackend.cpp:230
This class serves to eliminate the need to include the WellContributions into the matrix (with –matri...
Definition: WellContributions.hpp:53
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition: BlackoilPhases.hpp:27