cusparseSolverBackend.hpp
Go to the documentation of this file.
1/*
2 Copyright 2019 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_CUSPARSESOLVER_BACKEND_HEADER_INCLUDED
21#define OPM_CUSPARSESOLVER_BACKEND_HEADER_INCLUDED
22
23
24#include "cublas_v2.h"
25#include "cusparse_v2.h"
26
30
31namespace Opm
32{
33namespace Accelerator
34{
35
37template <unsigned int block_size>
38class cusparseSolverBackend : 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::deviceID;
48 using Base::maxit;
49 using Base::tolerance;
51
52private:
53
54 cublasHandle_t cublasHandle;
55 cusparseHandle_t cusparseHandle;
56 cudaStream_t stream;
57 cusparseMatDescr_t descr_B, descr_M, descr_L, descr_U;
58 bsrilu02Info_t info_M;
59 bsrsv2Info_t info_L, info_U;
60 // b: bsr matrix, m: preconditioner
61 double *d_bVals, *d_mVals;
62 int *d_bCols, *d_mCols;
63 int *d_bRows, *d_mRows;
64 double *d_x, *d_b, *d_r, *d_rw, *d_p; // vectors, used during linear solve
65 double *d_pw, *d_s, *d_t, *d_v;
66 void *d_buffer;
67 double *vals_contiguous; // only used if COPY_ROW_BY_ROW is true in cusparseSolverBackend.cpp
68
69 bool analysis_done = false;
70
71 bool useJacMatrix = false;
72 int nnzbs_prec; // number of nonzero blocks in the matrix for preconditioner
73 // could be jacMatrix or matrix
74
78 void gpu_pbicgstab(WellContributions& wellContribs, BdaResult& res);
79
83 void initialize(std::shared_ptr<BlockedMatrix> matrix, std::shared_ptr<BlockedMatrix> jacMatrix);
84
86 void finalize();
87
93 void copy_system_to_gpu(std::shared_ptr<BlockedMatrix> matrix, double *b, std::shared_ptr<BlockedMatrix> jacMatrix);
94
100 void update_system_on_gpu(std::shared_ptr<BlockedMatrix> matrix, double *b, std::shared_ptr<BlockedMatrix> jacMatrix);
101
104 bool analyse_matrix();
105
108 bool create_preconditioner();
109
113 void solve_system(WellContributions& wellContribs, BdaResult &res);
114
115public:
116
117
123 cusparseSolverBackend(int linear_solver_verbosity, int maxit, double tolerance, unsigned int deviceID);
124
127
135 SolverStatus solve_system(std::shared_ptr<BlockedMatrix> matrix, double *b,
136 std::shared_ptr<BlockedMatrix> jacMatrix, WellContributions& wellContribs, BdaResult &res) override;
137
140 void get_result(double *x) override;
141
142}; // end class cusparseSolverBackend
143
144} // namespace Accelerator
145} // namespace Opm
146
147#endif
148
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
This class implements a cusparse-based ilu0-bicgstab solver on GPU.
Definition: cusparseSolverBackend.hpp:38
~cusparseSolverBackend()
Destroy a cusparseSolver, and free memory.
SolverStatus solve_system(std::shared_ptr< BlockedMatrix > matrix, double *b, std::shared_ptr< BlockedMatrix > jacMatrix, WellContributions &wellContribs, BdaResult &res) override
cusparseSolverBackend(int linear_solver_verbosity, int maxit, double tolerance, unsigned int deviceID)
void get_result(double *x) override
Definition: WellContributions.hpp:52
SolverStatus
Definition: BdaSolver.hpp:35
Definition: BlackoilPhases.hpp:27