openclCPR.hpp
Go to the documentation of this file.
1/*
2 Copyright 2021 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_OPENCLCPR_HPP
21#define OPM_OPENCLCPR_HPP
22
23#include <dune/istl/paamg/matrixhierarchy.hh>
24#include <dune/istl/umfpack.hh>
25
32
34
35#include <type_traits>
36
37namespace Opm::Accelerator {
38
39template<class Scalar> class BlockedMatrix;
40
42template<class Scalar, unsigned int block_size>
43class openclCPR : public openclPreconditioner<Scalar,block_size>, public CprCreation<Scalar, block_size>
44{
46
47 using Base::N;
48 using Base::Nb;
49 using Base::nnz;
50 using Base::nnzb;
51 using Base::verbosity;
52 using Base::context;
53 using Base::queue;
54 using Base::events;
55 using Base::err;
56
57private:
58 std::vector<OpenclMatrix<Scalar>> d_Amatrices, d_Rmatrices; // scalar matrices that represent the AMG hierarchy
59 std::vector<cl::Buffer> d_PcolIndices;
60 std::vector<cl::Buffer> d_invDiags;
61 std::vector<cl::Buffer> d_t, d_f, d_u; // intermediate vectors used during amg cycle
62 std::unique_ptr<cl::Buffer> d_rs; // use before extracting the pressure
63 std::unique_ptr<cl::Buffer> d_weights; // the quasiimpes weights, used to extract pressure
64 std::unique_ptr<OpenclMatrix<Scalar>> d_mat; // stores blocked matrix
65 std::unique_ptr<cl::Buffer> d_coarse_y, d_coarse_x; // stores the scalar vectors
66 std::once_flag opencl_buffers_allocated; // only allocate OpenCL Buffers once
67
68 std::unique_ptr<openclBILU0<Scalar,block_size>> bilu0; // Blocked ILU0 preconditioner
69
70 std::unique_ptr<openclSolverBackend<Scalar,1>> coarse_solver; // coarse solver is scalar
71 bool opencl_ilu_parallel; // whether ILU0 operation should be parallelized
72
73 // Initialize and allocate matrices and vectors
74 void init_opencl_buffers();
75
76 // Copy matrices and vectors to GPU
77 void opencl_upload();
78
79 // apply pressure correction to vector
80 void apply_amg(const cl::Buffer& y, cl::Buffer& x);
81
82 void amg_cycle_gpu(const int level, cl::Buffer &y, cl::Buffer &x);
83
84public:
85 openclCPR(bool opencl_ilu_parallel, int verbosity);
86
89 BlockedMatrix<Scalar>* jacMat) override;
90
91 // set own Opencl variables, but also that of the bilu0 preconditioner
92 void setOpencl(std::shared_ptr<cl::Context>& context,
93 std::shared_ptr<cl::CommandQueue>& queue) override;
94
95 // applies blocked ilu0
96 // also applies amg for pressure component
97 void apply(const cl::Buffer& y, cl::Buffer& x) override;
98 void apply(Scalar&, Scalar&) override {}
99
102 BlockedMatrix<Scalar>* jacMat) override;
103};
104
105} // namespace Opm::Accelerator
106
107#endif
108
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 nnz
Definition: Preconditioner.hpp:45
int Nb
Definition: Preconditioner.hpp:44
int verbosity
Definition: Preconditioner.hpp:47
int nnzb
Definition: Preconditioner.hpp:46
int N
Definition: Preconditioner.hpp:43
This class implements a Constrained Pressure Residual (CPR) preconditioner.
Definition: openclCPR.hpp:44
void apply(const cl::Buffer &y, cl::Buffer &x) override
bool create_preconditioner(BlockedMatrix< Scalar > *mat, BlockedMatrix< Scalar > *jacMat) override
openclCPR(bool opencl_ilu_parallel, int verbosity)
bool analyze_matrix(BlockedMatrix< Scalar > *mat) override
bool analyze_matrix(BlockedMatrix< Scalar > *mat, BlockedMatrix< Scalar > *jacMat) override
bool create_preconditioner(BlockedMatrix< Scalar > *mat) override
void apply(Scalar &, Scalar &) override
Definition: openclCPR.hpp:98
void setOpencl(std::shared_ptr< cl::Context > &context, std::shared_ptr< cl::CommandQueue > &queue) override
Definition: openclPreconditioner.hpp:32
cl_int err
Definition: openclPreconditioner.hpp:38
std::shared_ptr< cl::Context > context
Definition: openclPreconditioner.hpp:35
std::vector< cl::Event > events
Definition: openclPreconditioner.hpp:37
std::shared_ptr< cl::CommandQueue > queue
Definition: openclPreconditioner.hpp:36
Definition: amgclSolverBackend.hpp:44