openclBISAI.hpp
Go to the documentation of this file.
1/*
2 Copyright 2022 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_OPENCLBISAI_HPP
21#define OPM_OPENCLBISAI_HPP
22
23#include <mutex>
24
28
29namespace Opm::Accelerator {
30
31template<class Scalar> class BlockedMatrix;
32
35template<class Scalar, unsigned int block_size>
36class openclBISAI : public openclPreconditioner<Scalar,block_size>
37{
39
40 using Base::N;
41 using Base::Nb;
42 using Base::nnz;
43 using Base::nnzb;
44 using Base::verbosity;
45 using Base::context;
46 using Base::queue;
47 using Base::events;
48 using Base::err;
49
50private:
51 std::once_flag initialize;
52
53 std::vector<int> colPointers;
54 std::vector<int> rowIndices;
55 std::vector<int> diagIndex;
56 std::vector<int> csrToCscOffsetMap;
57 std::vector<Scalar> invLvals;
58 std::vector<Scalar> invUvals;
59
60 cl::Buffer d_colPointers;
61 cl::Buffer d_rowIndices;
62 cl::Buffer d_csrToCscOffsetMap;
63 cl::Buffer d_diagIndex;
64 cl::Buffer d_LUvals;
65 cl::Buffer d_invDiagVals;
66 cl::Buffer d_invLvals;
67 cl::Buffer d_invUvals;
68 cl::Buffer d_invL_x;
69
70 bool opencl_ilu_parallel;
71 std::unique_ptr<openclBILU0<Scalar,block_size>> bilu0;
72
74 struct subsystemStructure {
77 std::vector<int> subsystemPointers;
82 std::vector<int> nzIndices;
85 std::vector<int> knownRhsIndices;
87 std::vector<int> unknownRhsIndices;
88 };
89
91 struct subsystemStructureGPU {
92 cl::Buffer subsystemPointers;
93 cl::Buffer nzIndices;
94 cl::Buffer knownRhsIndices;
95 cl::Buffer unknownRhsIndices;
96 } ;
97
98 subsystemStructure lower, upper;
99 subsystemStructureGPU d_lower, d_upper;
100
103 void buildLowerSubsystemsStructures();
104
107 void buildUpperSubsystemsStructures();
108
109public:
110 openclBISAI(bool opencl_ilu_parallel, int verbosity);
111
112 // set own Opencl variables, but also that of the bilu0 preconditioner
113 void setOpencl(std::shared_ptr<cl::Context>& context,
114 std::shared_ptr<cl::CommandQueue>& queue) override;
115
116 // analysis, extract parallelism
119 BlockedMatrix<Scalar>* jacMat) override;
120
121 // ilu_decomposition
124 BlockedMatrix<Scalar>* jacMat) override;
125
126 // apply preconditioner, x = prec(y)
127 void apply(const cl::Buffer& y, cl::Buffer& x) override;
128 void apply(Scalar&, Scalar&) override {}
129};
130
134std::vector<int> buildCsrToCscOffsetMap(std::vector<int> colPointers, std::vector<int> rowIndices);
135
136} // namespace Opm::Accelerator
137
138#endif
Definition: BlockedMatrix.hpp:29
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
Definition: openclBISAI.hpp:37
bool create_preconditioner(BlockedMatrix< Scalar > *mat) override
bool analyze_matrix(BlockedMatrix< Scalar > *mat, BlockedMatrix< Scalar > *jacMat) override
void apply(const cl::Buffer &y, cl::Buffer &x) override
void apply(Scalar &, Scalar &) override
Definition: openclBISAI.hpp:128
bool create_preconditioner(BlockedMatrix< Scalar > *mat, BlockedMatrix< Scalar > *jacMat) override
openclBISAI(bool opencl_ilu_parallel, int verbosity)
void setOpencl(std::shared_ptr< cl::Context > &context, std::shared_ptr< cl::CommandQueue > &queue) override
bool analyze_matrix(BlockedMatrix< Scalar > *mat) 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
std::vector< int > buildCsrToCscOffsetMap(std::vector< int > colPointers, std::vector< int > rowIndices)