GpuDILU.hpp
Go to the documentation of this file.
1/*
2 Copyright 2022-2023 SINTEF AS
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#ifndef OPM_GPUDILU_HPP
20#define OPM_GPUDILU_HPP
21
22#include <memory>
23#include <opm/grid/utility/SparseTable.hpp>
28#include <vector>
29#include <map>
30#include <utility>
31
32
33namespace Opm::gpuistl
34{
48template <class CPUMatrixT, class X, class Y, int l = 1>
50{
51public:
53 using domain_type = X;
55 using range_type = Y;
57 using field_type = typename X::field_type;
58
61 using FloatVec = GpuVector<float>;
62
65
66
70 explicit GpuDILU(const GPUMatrix& gpuMatrix, const CPUMatrixT& cpuMatrix, bool splitMatrix, bool tuneKernels, int mixedPrecisionScheme);
71
74 void pre(X& x, Y& b) override;
75
77 void apply(X& v, const Y& d) override;
78
81 void post(X& x) override;
82
84 Dune::SolverCategory::Category category() const override;
85
87 void update() final;
88
90 void reorderAndSplitMatrix(int moveThreadBlockSize);
91
93 void computeDiagonal(int factorizationThreadBlockSize);
94
97
98
100 static constexpr bool shouldCallPre()
101 {
102 return false;
103 }
104
106 static constexpr bool shouldCallPost()
107 {
108 return false;
109 }
110
111 virtual bool hasPerfectUpdate() const override {
112 return true;
113 }
114
115
116private:
118 void apply(X& v, const Y& d, int lowerSolveThreadBlockSize, int upperSolveThreadBlockSize);
120 void update(int moveThreadBlockSize, int factorizationThreadBlockSize);
122 static constexpr const size_t blocksize_ = CPUMatrixT::block_type::cols;
124 Opm::SparseTable<size_t> m_levelSets;
126 std::vector<int> m_reorderedToNatural;
128 std::vector<int> m_naturalToReordered;
130 const GPUMatrix& m_gpuMatrix;
132 std::unique_ptr<GPUMatrix> m_gpuMatrixReordered;
134 std::unique_ptr<GPUMatrix> m_gpuMatrixReorderedLower;
135 std::unique_ptr<GPUMatrix> m_gpuMatrixReorderedUpper;
137 std::unique_ptr<GpuVector<field_type>> m_gpuMatrixReorderedDiag;
139 std::unique_ptr<FloatMat> m_gpuMatrixReorderedLowerFloat;
140 std::unique_ptr<FloatMat> m_gpuMatrixReorderedUpperFloat;
141 std::unique_ptr<FloatVec> m_gpuMatrixReorderedDiagFloat;
142 std::unique_ptr<FloatVec> m_gpuDInvFloat;
144 GpuVector<int> m_gpuNaturalToReorder;
146 GpuVector<int> m_gpuReorderToNatural;
148 GpuVector<field_type> m_gpuDInv;
150 bool m_splitMatrix;
152 bool m_tuneThreadBlockSizes;
154 const MatrixStorageMPScheme m_mixedPrecisionScheme;
157 int m_upperSolveThreadBlockSize = -1;
158 int m_lowerSolveThreadBlockSize = -1;
159 int m_moveThreadBlockSize = -1;
160 int m_DILUFactorizationThreadBlockSize = -1;
161
162 // Graphs for Apply
163 std::map<std::pair<field_type*, const field_type*>, GPUGraph> m_apply_graphs;
164 std::map<std::pair<field_type*, const field_type*>, GPUGraphExec> m_executableGraphs;
165
166 // Stream for the DILU operations on the GPU
167 GPUStream m_stream{};
168 // Events for synchronization with main stream
169 GPUEvent m_before{};
170 GPUEvent m_after{};
171};
172} // end namespace Opm::gpuistl
173
174#endif
Interface class adding the update() method to the preconditioner interface.
Definition: PreconditionerWithUpdate.hpp:32
DILU preconditioner on the GPU.
Definition: GpuDILU.hpp:50
void pre(X &x, Y &b) override
Prepare the preconditioner.
void tuneThreadBlockSizes()
function that will experimentally tune the thread block sizes of the important cuda kernels
void update() final
Updates the matrix data.
void apply(X &v, const Y &d) override
Apply the preconditoner.
GpuDILU(const GPUMatrix &gpuMatrix, const CPUMatrixT &cpuMatrix, bool splitMatrix, bool tuneKernels, int mixedPrecisionScheme)
Constructor.
virtual bool hasPerfectUpdate() const override
Definition: GpuDILU.hpp:111
Y range_type
The range type of the preconditioner.
Definition: GpuDILU.hpp:55
static constexpr bool shouldCallPre()
Definition: GpuDILU.hpp:100
Dune::SolverCategory::Category category() const override
Category of the preconditioner (see SolverCategory::Category)
GpuVector< float > FloatVec
Definition: GpuDILU.hpp:61
void computeDiagonal(int factorizationThreadBlockSize)
Compute the diagonal of the DILU, and update the data of the reordered matrix.
static constexpr bool shouldCallPost()
Definition: GpuDILU.hpp:106
X domain_type
The domain type of the preconditioner.
Definition: GpuDILU.hpp:53
void reorderAndSplitMatrix(int moveThreadBlockSize)
perform matrix splitting and reordering
typename X::field_type field_type
The field type of the preconditioner.
Definition: GpuDILU.hpp:57
void post(X &x) override
Post processing.
GpuSparseMatrix< field_type > GPUMatrix
Definition: GpuDILU.hpp:59
Definition: autotuner.hpp:29
MatrixStorageMPScheme
Definition: kernel_enums.hpp:31