My Project
Loading...
Searching...
No Matches
rocsparseCPR.hpp
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
25#include <opm/simulators/linalg/bda/rocm/rocsparseBILU0.hpp>
26#include <opm/simulators/linalg/bda/Matrix.hpp>
27#include <opm/simulators/linalg/bda/CprCreation.hpp>
28#include <opm/simulators/linalg/bda/rocm/rocsparseMatrix.hpp>
29#include <opm/simulators/linalg/bda/rocm/rocsparsePreconditioner.hpp>
30
31#include <opm/simulators/linalg/bda/rocm/rocsparseSolverBackend.hpp>
32
33namespace Opm::Accelerator {
34
35template<class Scalar> class BlockedMatrix;
36
38template <class Scalar, unsigned int block_size>
39class rocsparseCPR : public rocsparsePreconditioner<Scalar, block_size>, public CprCreation<Scalar, block_size>
40{
42
43 using Base::N;
44 using Base::Nb;
45 using Base::nnz;
46 using Base::nnzb;
47 using Base::verbosity;
48
49private:
50 std::vector<RocmMatrix<Scalar>> d_Amatrices, d_Rmatrices; // scalar matrices that represent the AMG hierarchy
51
52 std::vector<RocmVector<int>> d_PcolIndices; // prolongation does not need a full matrix, only store colIndices
53 std::vector<RocmVector<Scalar>> d_invDiags; // inverse of diagonal of Amatrices
54 std::vector<RocmVector<Scalar>> d_t, d_f; // intermediate vectors used during amg cycle
55 std::vector<RocmVector<Scalar>> d_u; // intermediate vectors used during amg cycle
56 std::vector<RocmVector<Scalar>> d_rs; // use before extracting the pressure
57 std::vector<RocmVector<Scalar>> d_weights; // the quasiimpes weights, used to extract pressure
58 std::unique_ptr<RocmMatrix<Scalar>> d_mat; // stores blocked matrix
59 std::vector<RocmVector<Scalar>> d_coarse_y, d_coarse_x; // stores the scalar vectors
60 std::once_flag rocm_buffers_allocated; // only allocate OpenCL Buffers once
61
62 std::unique_ptr<rocsparseBILU0<Scalar, block_size> > bilu0; // Blocked ILU0 preconditioner
63
64 std::unique_ptr<rocsparseSolverBackend<Scalar, 1> > coarse_solver; // coarse solver is scalar
65
66 // Initialize and allocate matrices and vectors
67 void init_rocm_buffers();
68
69 // Copy matrices and vectors to GPU
70 void rocm_upload();
71
72 // apply pressure correction to vector
73 void apply_amg(const Scalar& y, Scalar& x);
74
75 // Apply the AMG preconditioner
76 void amg_cycle_gpu(const int level, Scalar &y, Scalar &x);
77
78public:
79
80 rocsparseCPR(int verbosity);
81
85 bool initialize(std::shared_ptr<BlockedMatrix<Scalar>> matrix,
86 std::shared_ptr<BlockedMatrix<Scalar>> jacMatrix,
87 rocsparse_int *d_Arows,
88 rocsparse_int *d_Acols) override;
89
90
93 bool analyze_matrix(BlockedMatrix<Scalar> *mat) override;
94
99 BlockedMatrix<Scalar> *jacMat) override;
100
104
109 BlockedMatrix<Scalar> *jacMat) override;
110
116 void apply(Scalar& y,
117 Scalar& x) override;
118
119#if HAVE_OPENCL
120 // apply preconditioner, x = prec(y)
121 void apply(const cl::Buffer&, cl::Buffer&) override {}
122#endif
123
126 void copy_system_to_gpu(Scalar *b) override;
127
131// void update_system(Scalar *vals, Scalar *b);
132
135 void update_system_on_gpu(Scalar *b) override;
136};
137
138} // namespace Opm
139
140#endif
141
This struct resembles a blocked csr matrix, like Dune::BCRSMatrix.
Definition BlockedMatrix.hpp:29
This class implements a Constrained Pressure Residual (CPR) preconditioner.
Definition CprCreation.hpp:39
This class implements a Constrained Pressure Residual (CPR) preconditioner.
Definition rocsparseCPR.hpp:40
bool analyze_matrix(BlockedMatrix< Scalar > *mat) override
Analysis, extract parallelism if specified.
Definition rocsparseCPR.cpp:84
bool initialize(std::shared_ptr< BlockedMatrix< Scalar > > matrix, std::shared_ptr< BlockedMatrix< Scalar > > jacMatrix, rocsparse_int *d_Arows, rocsparse_int *d_Acols) override
Initialize GPU and allocate memory.
Definition rocsparseCPR.cpp:54
void apply(Scalar &y, Scalar &x) override
Apply preconditioner, x = prec(y) applies blocked ilu0 also applies amg for pressure component.
Definition rocsparseCPR.cpp:309
void copy_system_to_gpu(Scalar *b) override
Copy matrix A values to GPU.
Definition rocsparseCPR.cpp:72
bool create_preconditioner(BlockedMatrix< Scalar > *mat) override
Create AMG preconditioner and perform ILU decomposition.
Definition rocsparseCPR.cpp:146
void update_system_on_gpu(Scalar *b) override
Reassign pointers, in case the addresses of the Dune variables have changed --> TODO: check when/if w...
Definition rocsparseCPR.cpp:78
Definition rocsparsePreconditioner.hpp:33
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:242