My Project
Loading...
Searching...
No Matches
GpuDILU.hpp
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>
24#include <opm/simulators/linalg/PreconditionerWithUpdate.hpp>
25#include <opm/simulators/linalg/gpuistl/GpuSparseMatrix.hpp>
26#include <opm/simulators/linalg/gpuistl/detail/kernel_enums.hpp>
27#include <opm/simulators/linalg/gpuistl/gpu_resources.hpp>
28#include <vector>
29#include <map>
30#include <utility>
31
32
33namespace Opm::gpuistl
34{
45template <class M, class X, class Y, int l = 1>
47{
48public:
50 using matrix_type = typename std::remove_const<M>::type;
52 using domain_type = X;
54 using range_type = Y;
56 using field_type = typename X::field_type;
60 using FloatVec = GpuVector<float>;
61
68 explicit GpuDILU(const M& A, bool splitMatrix, bool tuneKernels, int mixedPrecisionScheme);
69
72 void pre(X& x, Y& b) override;
73
75 void apply(X& v, const Y& d) override;
76
79 void post(X& x) override;
80
82 Dune::SolverCategory::Category category() const override;
83
85 void update() final;
86
89
92
95
96
99 {
100 return false;
101 }
102
104 static constexpr bool shouldCallPost()
105 {
106 return false;
107 }
108
109 virtual bool hasPerfectUpdate() const override {
110 return true;
111 }
112
113
114private:
116 void apply(X& v, const Y& d, int lowerSolveThreadBlockSize, int upperSolveThreadBlockSize);
120 const M& m_cpuMatrix;
122 static constexpr const size_t blocksize_ = matrix_type::block_type::cols;
124 Opm::SparseTable<size_t> m_levelSets;
126 std::vector<int> m_reorderedToNatural;
128 std::vector<int> m_naturalToReordered;
130 CuMat m_gpuMatrix;
132 std::unique_ptr<CuMat> m_gpuMatrixReordered;
134 std::unique_ptr<CuMat> m_gpuMatrixReorderedLower;
135 std::unique_ptr<CuMat> 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
Definition BlackoilWellModel.hpp:83
DILU preconditioner on the GPU.
Definition GpuDILU.hpp:47
void apply(X &v, const Y &d) override
Apply the preconditoner.
Definition GpuDILU.cpp:113
Y range_type
The range type of the preconditioner.
Definition GpuDILU.hpp:54
void computeDiagonal(int factorizationThreadBlockSize)
Compute the diagonal of the DILU, and update the data of the reordered matrix.
Definition GpuDILU.cpp:339
GpuSparseMatrix< field_type > CuMat
The GPU matrix type.
Definition GpuDILU.hpp:58
typename std::remove_const< M >::type matrix_type
The matrix type the preconditioner is for.
Definition GpuDILU.hpp:50
void pre(X &x, Y &b) override
Prepare the preconditioner.
Definition GpuDILU.cpp:107
void post(X &x) override
Post processing.
Definition GpuDILU.cpp:269
typename X::field_type field_type
The field type of the preconditioner.
Definition GpuDILU.hpp:56
static constexpr bool shouldCallPre()
Definition GpuDILU.hpp:98
Dune::SolverCategory::Category category() const override
Category of the preconditioner (see SolverCategory::Category)
Definition GpuDILU.cpp:275
void reorderAndSplitMatrix(int moveThreadBlockSize)
perform matrix splitting and reordering
Definition GpuDILU.cpp:311
void tuneThreadBlockSizes()
function that will experimentally tune the thread block sizes of the important cuda kernels
Definition GpuDILU.cpp:419
X domain_type
The domain type of the preconditioner.
Definition GpuDILU.hpp:52
static constexpr bool shouldCallPost()
Definition GpuDILU.hpp:104
void update() final
Updates the matrix data.
Definition GpuDILU.cpp:282
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:242