My Project
Loading...
Searching...
No Matches
rocsparseSolverBackend.hpp
1/*
2 Copyright 2023 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_ROCSPARSESOLVER_BACKEND_HEADER_INCLUDED
21#define OPM_ROCSPARSESOLVER_BACKEND_HEADER_INCLUDED
22
23#include <memory>
24
25#include <opm/simulators/linalg/gpubridge/GpuResult.hpp>
26#include <opm/simulators/linalg/gpubridge/GpuSolver.hpp>
27#include <opm/simulators/linalg/gpubridge/WellContributions.hpp>
28
29#include <opm/simulators/linalg/gpubridge/rocm/rocsparsePreconditioner.hpp>
30
31#include <rocblas/rocblas.h>
32#include <rocsparse/rocsparse.h>
33
34#include <hip/hip_version.h>
35
36namespace Opm::Accelerator {
37
39template<class Scalar, unsigned int block_size>
40class rocsparseSolverBackend : public GpuSolver<Scalar,block_size>
41{
43
44 using Base::N;
45 using Base::Nb;
46 using Base::nnz;
47 using Base::nnzb;
48 using Base::verbosity;
49 using Base::platformID;
50 using Base::deviceID;
51 using Base::maxit;
52 using Base::tolerance;
53 using Base::initialized;
54
55private:
56 double c_copy = 0.0; // cummulative timer measuring the total time it takes to transfer the data to the GPU
57
58 bool useJacMatrix = false;
59
60 bool analysis_done = false;
61 std::shared_ptr<BlockedMatrix<Scalar>> mat{}; // original matrix
62
65 rocsparse_handle handle;
66 rocblas_handle blas_handle;
67 rocsparse_mat_descr descr_A;
68#if HIP_VERSION >= 50400000
70#endif
71 hipStream_t stream;
72
73 rocsparse_int *d_Arows, *d_Acols;
74 Scalar *d_Avals;
75 Scalar *d_x, *d_b, *d_r, *d_rw, *d_p; // vectors, used during linear solve
76 Scalar *d_pw, *d_s, *d_t, *d_v;
77 int ver;
78 char rev[64];
79
80 std::unique_ptr<rocsparsePreconditioner<Scalar, block_size> > prec; // can perform blocked ILU0 and AMG on pressure component
81
86
90 void initialize(std::shared_ptr<BlockedMatrix<Scalar>> matrix,
91 std::shared_ptr<BlockedMatrix<Scalar>> jacMatrix);
92
95 void copy_system_to_gpu(Scalar* b);
96
99 void update_system_on_gpu(Scalar* b);
100
103 bool analyze_matrix();
104
107 bool create_preconditioner();
108
113
114public:
123 int maxit,
124 Scalar tolerance,
125 unsigned int platformID,
126 unsigned int deviceID,
127 std::string linsolver);
128
131
134
142 SolverStatus solve_system(std::shared_ptr<BlockedMatrix<Scalar>> matrix,
143 Scalar* b,
144 std::shared_ptr<BlockedMatrix<Scalar>> jacMatrix,
146 GpuResult& res) override;
147
150 void get_result(Scalar* x) override;
151
152}; // end class rocsparseSolverBackend
153
154} // namespace Opm::Accelerator
155
156#endif
This struct resembles a blocked csr matrix, like Dune::BCRSMatrix.
Definition BlockedMatrix.hpp:29
This class is based on InverseOperatorResult struct from dune/istl/solver.hh It is needed to prevent ...
Definition GpuResult.hpp:31
This class serves to simplify choosing between different backend solvers, such as cusparseSolver and ...
Definition GpuSolver.hpp:46
This class implements a rocsparse-based ilu0-bicgstab solver on GPU.
Definition rocsparseSolverBackend.hpp:41
~rocsparseSolverBackend()
Destroy a openclSolver, and free memory.
Definition rocsparseSolverBackend.cpp:113
void get_result(Scalar *x) override
Get result after linear solve, and peform postprocessing if necessary.
Definition rocsparseSolverBackend.cpp:660
rocsparseSolverBackend(int linear_solver_verbosity, int maxit, Scalar tolerance, bool opencl_ilu_reorder)
For the CPR coarse solver.
This class serves to eliminate the need to include the WellContributions into the matrix (with –matri...
Definition WellContributions.hpp:51
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:242