BISAI.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 BISAI_HPP
21#define BISAI_HPP
22
23#include <mutex>
24
28
29namespace Opm
30{
31namespace Accelerator
32{
33
34class BlockedMatrix;
35
38template <unsigned int block_size>
39class BISAI : public Preconditioner<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 using Base::context;
49 using Base::queue;
50 using Base::events;
51 using Base::err;
52
53private:
54 std::once_flag initialize;
55
56 std::vector<int> colPointers;
57 std::vector<int> rowIndices;
58 std::vector<int> diagIndex;
59 std::vector<int> csrToCscOffsetMap;
60 std::vector<double> invLvals;
61 std::vector<double> invUvals;
62
63 cl::Buffer d_colPointers;
64 cl::Buffer d_rowIndices;
65 cl::Buffer d_csrToCscOffsetMap;
66 cl::Buffer d_diagIndex;
67 cl::Buffer d_LUvals;
68 cl::Buffer d_invDiagVals;
69 cl::Buffer d_invLvals;
70 cl::Buffer d_invUvals;
71 cl::Buffer d_invL_x;
72
73 bool opencl_ilu_parallel;
74 std::unique_ptr<BILU0<block_size> > bilu0;
75
77 typedef struct{
80 std::vector<int> subsystemPointers;
85 std::vector<int> nzIndices;
88 std::vector<int> knownRhsIndices;
90 std::vector<int> unknownRhsIndices;
91 } subsystemStructure;
92
94 typedef struct{
95 cl::Buffer subsystemPointers;
96 cl::Buffer nzIndices;
97 cl::Buffer knownRhsIndices;
98 cl::Buffer unknownRhsIndices;
99 } subsystemStructureGPU;
100
101 subsystemStructure lower, upper;
102 subsystemStructureGPU d_lower, d_upper;
103
106 void buildLowerSubsystemsStructures();
107
110 void buildUpperSubsystemsStructures();
111
112public:
113 BISAI(bool opencl_ilu_parallel, int verbosity);
114
115 // set own Opencl variables, but also that of the bilu0 preconditioner
116 void setOpencl(std::shared_ptr<cl::Context>& context, std::shared_ptr<cl::CommandQueue>& queue) override;
117
118 // analysis, extract parallelism
119 bool analyze_matrix(BlockedMatrix *mat) override;
120 bool analyze_matrix(BlockedMatrix *mat, BlockedMatrix *jacMat) override;
121
122 // ilu_decomposition
125
126 // apply preconditioner, x = prec(y)
127 void apply(const cl::Buffer& y, cl::Buffer& x) override;
128};
129
133std::vector<int> buildCsrToCscOffsetMap(std::vector<int> colPointers, std::vector<int> rowIndices);
134
135} // namespace Accelerator
136} // namespace Opm
137
138#endif
Definition: BISAI.hpp:40
bool create_preconditioner(BlockedMatrix *mat) override
void apply(const cl::Buffer &y, cl::Buffer &x) override
void setOpencl(std::shared_ptr< cl::Context > &context, std::shared_ptr< cl::CommandQueue > &queue) override
BISAI(bool opencl_ilu_parallel, int verbosity)
bool create_preconditioner(BlockedMatrix *mat, BlockedMatrix *jacMat) override
bool analyze_matrix(BlockedMatrix *mat, BlockedMatrix *jacMat) override
bool analyze_matrix(BlockedMatrix *mat) override
Definition: BlockedMatrix.hpp:31
Definition: Preconditioner.hpp:36
int nnz
Definition: Preconditioner.hpp:41
std::vector< cl::Event > events
Definition: Preconditioner.hpp:47
std::shared_ptr< cl::CommandQueue > queue
Definition: Preconditioner.hpp:46
int N
Definition: Preconditioner.hpp:39
int verbosity
Definition: Preconditioner.hpp:43
cl_int err
Definition: Preconditioner.hpp:48
int Nb
Definition: Preconditioner.hpp:40
std::shared_ptr< cl::Context > context
Definition: Preconditioner.hpp:45
int nnzb
Definition: Preconditioner.hpp:42
std::vector< int > buildCsrToCscOffsetMap(std::vector< int > colPointers, std::vector< int > rowIndices)
Definition: BlackoilPhases.hpp:27