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::Accelerator {
32
34template<class Scalar, unsigned int block_size>
35class cusparseSolverBackend : public BdaSolver<Scalar,block_size>
36{
38
39 using Base::N;
40 using Base::Nb;
41 using Base::nnz;
42 using Base::nnzb;
43 using Base::verbosity;
44 using Base::deviceID;
45 using Base::maxit;
46 using Base::tolerance;
48
49private:
50 cublasHandle_t cublasHandle;
51 cusparseHandle_t cusparseHandle;
52 cudaStream_t stream;
53 cusparseMatDescr_t descr_B, descr_M, descr_L, descr_U;
54 bsrilu02Info_t info_M;
55 bsrsv2Info_t info_L, info_U;
56 // b: bsr matrix, m: preconditioner
57 Scalar *d_bVals, *d_mVals;
58 int *d_bCols, *d_mCols;
59 int *d_bRows, *d_mRows;
60 Scalar *d_x, *d_b, *d_r, *d_rw, *d_p; // vectors, used during linear solve
61 Scalar *d_pw, *d_s, *d_t, *d_v;
62 void *d_buffer;
63 Scalar *vals_contiguous; // only used if COPY_ROW_BY_ROW is true in cusparseSolverBackend.cpp
64
65 bool analysis_done = false;
66
67 bool useJacMatrix = false;
68 int nnzbs_prec; // number of nonzero blocks in the matrix for preconditioner
69 // could be jacMatrix or matrix
70
71 double c_copy = 0.0; // cummulative timer measuring the total time it takes to transfer the data to the GPU
72
76 void gpu_pbicgstab(WellContributions<Scalar>& wellContribs, BdaResult& res);
77
81 void initialize(std::shared_ptr<BlockedMatrix<Scalar>> matrix,
82 std::shared_ptr<BlockedMatrix<Scalar>> jacMatrix);
83
85 void finalize();
86
92 void copy_system_to_gpu(std::shared_ptr<BlockedMatrix<Scalar>> matrix,
93 Scalar* b,
94 std::shared_ptr<BlockedMatrix<Scalar>> jacMatrix);
95
101 void update_system_on_gpu(std::shared_ptr<BlockedMatrix<Scalar>> matrix,
102 Scalar* b,
103 std::shared_ptr<BlockedMatrix<Scalar>> jacMatrix);
104
107 bool analyse_matrix();
108
111 bool create_preconditioner();
112
116 void solve_system(WellContributions<Scalar>& wellContribs, BdaResult &res);
117
118public:
124 cusparseSolverBackend(int linear_solver_verbosity, int maxit,
125 Scalar tolerance, unsigned int deviceID);
126
129
138 Scalar* b,
139 std::shared_ptr<BlockedMatrix<Scalar>> jacMatrix,
140 WellContributions<Scalar>& wellContribs,
141 BdaResult& res) override;
142
145 void get_result(Scalar* x) override;
146
147}; // end class cusparseSolverBackend
148
149} // namespace Opm::Accelerator
150
151#endif
152
Definition: BdaResult.hpp:31
Definition: BdaSolver.hpp:46
int maxit
Definition: BdaSolver.hpp:55
int nnzb
Definition: BdaSolver.hpp:61
bool initialized
Definition: BdaSolver.hpp:66
int verbosity
Definition: BdaSolver.hpp:53
unsigned int deviceID
Definition: BdaSolver.hpp:64
int N
Definition: BdaSolver.hpp:58
int Nb
Definition: BdaSolver.hpp:59
int nnz
Definition: BdaSolver.hpp:60
Scalar tolerance
Definition: BdaSolver.hpp:56
Definition: BlockedMatrix.hpp:29
This class implements a cusparse-based ilu0-bicgstab solver on GPU.
Definition: cusparseSolverBackend.hpp:36
cusparseSolverBackend(int linear_solver_verbosity, int maxit, Scalar tolerance, unsigned int deviceID)
SolverStatus solve_system(std::shared_ptr< BlockedMatrix< Scalar > > matrix, Scalar *b, std::shared_ptr< BlockedMatrix< Scalar > > jacMatrix, WellContributions< Scalar > &wellContribs, BdaResult &res) override
~cusparseSolverBackend()
Destroy a cusparseSolver, and free memory.
void get_result(Scalar *x) override
Definition: WellContributions.hpp:53
Definition: amgclSolverBackend.hpp:44
SolverStatus
Definition: BdaSolver.hpp:35