rocsparseCPR.hpp
Go to the documentation of this file.
1/*
2 Copyright 2024 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_ROCSPARSECPR_HPP
21#define OPM_ROCSPARSECPR_HPP
22
23#include <mutex>
24
31
33
34namespace Opm::Accelerator {
35
36template<class Scalar> class BlockedMatrix;
37
39template <class Scalar, unsigned int block_size>
40class rocsparseCPR : public rocsparsePreconditioner<Scalar, block_size>, public CprCreation<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
50private:
51 std::vector<RocmMatrix<Scalar>> d_Amatrices, d_Rmatrices; // scalar matrices that represent the AMG hierarchy
52
53 std::vector<RocmVector<int>> d_PcolIndices; // prolongation does not need a full matrix, only store colIndices
54 std::vector<RocmVector<Scalar>> d_invDiags; // inverse of diagonal of Amatrices
55 std::vector<RocmVector<Scalar>> d_t, d_f; // intermediate vectors used during amg cycle
56 std::vector<RocmVector<Scalar>> d_u; // intermediate vectors used during amg cycle
57 std::vector<RocmVector<Scalar>> d_rs; // use before extracting the pressure
58 std::vector<RocmVector<Scalar>> d_weights; // the quasiimpes weights, used to extract pressure
59 std::unique_ptr<RocmMatrix<Scalar>> d_mat; // stores blocked matrix
60 std::shared_ptr<BlockedMatrix<Scalar>> h_mat;
61
62 std::vector<RocmVector<Scalar>> d_coarse_y, d_coarse_x, d_yinput; // stores the scalar vectors
63 std::once_flag rocm_buffers_allocated; // only allocate OpenCL Buffers once
64
65 std::unique_ptr<rocsparseBILU0<Scalar, block_size> > bilu0; // Blocked ILU0 preconditioner
66
67 std::vector<std::shared_ptr<rocsparseBILU0<Scalar,1>>> smootherilu0s;
68
69 Scalar *d_tmp, *d_yamg, *d_ywell;
70 bool createfirsttime;
71
72 #if HIP_VERSION >= 50400000
73 rocsparse_mat_info spmv_info;
74 #endif
75
76 // Initialize and allocate matrices and vectors
77 void init_rocm_buffers();
78 void init_rocm_buffers_fineMatrix();
79
80 // Deallocate matrices and vectors
81 void clear_rocm_buffers();
82
83 // Copy matrices and vectors to GPU
84 void rocm_upload();
85
86 // apply pressure correction to vector
87 void apply_amg(const Scalar& y,
88 Scalar& x,
89 WellContributions<Scalar>& wellContribs);
90
91 // Apply the AMG preconditioner
92 void amg_cycle_gpu(const int level,
93 Scalar &y,
94 Scalar &x,
95 WellContributions<Scalar>& wellContribs);
96
97public:
98
100
106 bool initialize(std::shared_ptr<BlockedMatrix<Scalar>> matrix,
107 std::shared_ptr<BlockedMatrix<Scalar>> jacMatrix,
108 rocsparse_int *d_Arows,
109 rocsparse_int *d_Acols) override;
110
111
115
120 BlockedMatrix<Scalar> *jacMat) override;
121
125
130 BlockedMatrix<Scalar> *jacMat) override;
131
138 void apply(const Scalar& y,
139 Scalar& x,
140 WellContributions<Scalar>& wellContribs) override;
141
144 void copy_system_to_gpu(Scalar *b) override;
145
149 void update_system_on_gpu(Scalar* vals, Scalar* b) override;
150
151};
152
153} // namespace Opm
154
155#endif
156
Definition: BlockedMatrix.hpp:29
This class implements a Constrained Pressure Residual (CPR) preconditioner.
Definition: CprCreation.hpp:39
BlockedMatrix< Scalar > * mat
Definition: CprCreation.hpp:56
int nnzb
Definition: Preconditioner.hpp:48
int Nb
Definition: Preconditioner.hpp:46
int N
Definition: Preconditioner.hpp:45
int nnz
Definition: Preconditioner.hpp:47
int verbosity
Definition: Preconditioner.hpp:49
This class implements a Constrained Pressure Residual (CPR) preconditioner.
Definition: rocsparseCPR.hpp:41
void update_system_on_gpu(Scalar *vals, Scalar *b) override
bool analyze_matrix(BlockedMatrix< Scalar > *mat, BlockedMatrix< Scalar > *jacMat) override
bool analyze_matrix(BlockedMatrix< Scalar > *mat) override
bool initialize(std::shared_ptr< BlockedMatrix< Scalar > > matrix, std::shared_ptr< BlockedMatrix< Scalar > > jacMatrix, rocsparse_int *d_Arows, rocsparse_int *d_Acols) override
void apply(const Scalar &y, Scalar &x, WellContributions< Scalar > &wellContribs) override
bool create_preconditioner(BlockedMatrix< Scalar > *mat, BlockedMatrix< Scalar > *jacMat) override
void copy_system_to_gpu(Scalar *b) override
bool create_preconditioner(BlockedMatrix< Scalar > *mat) override
Definition: rocsparsePreconditioner.hpp:34
std::shared_ptr< BlockedMatrix< Scalar > > jacMat
Definition: rocsparsePreconditioner.hpp:53
Definition: WellContributions.hpp:51
Definition: amgclSolverBackend.hpp:44