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>
26#include <vector>
27
28
29
30namespace Opm::gpuistl
31{
42template <class M, class X, class Y, int l = 1>
44{
45public:
47 using matrix_type = typename std::remove_const<M>::type;
49 using domain_type = X;
51 using range_type = Y;
53 using field_type = typename X::field_type;
56
63 explicit GpuDILU(const M& A, bool splitMatrix, bool tuneKernels);
64
67 void pre(X& x, Y& b) override;
68
70 void apply(X& v, const Y& d) override;
71
74 void post(X& x) override;
75
77 Dune::SolverCategory::Category category() const override;
78
80 void update() final;
81
83 void computeDiagAndMoveReorderedData(int moveThreadBlockSize, int factorizationThreadBlockSize);
84
87
88
90 static constexpr bool shouldCallPre()
91 {
92 return false;
93 }
94
96 static constexpr bool shouldCallPost()
97 {
98 return false;
99 }
100
101 virtual bool hasPerfectUpdate() const override {
102 return true;
103 }
104
105
106private:
108 void apply(X& v, const Y& d, int lowerSolveThreadBlockSize, int upperSolveThreadBlockSize);
110 void update(int moveThreadBlockSize, int factorizationThreadBlockSize);
112 const M& m_cpuMatrix;
114 static constexpr const size_t blocksize_ = matrix_type::block_type::cols;
116 Opm::SparseTable<size_t> m_levelSets;
118 std::vector<int> m_reorderedToNatural;
120 std::vector<int> m_naturalToReordered;
122 CuMat m_gpuMatrix;
124 std::unique_ptr<CuMat> m_gpuMatrixReordered;
126 std::unique_ptr<CuMat> m_gpuMatrixReorderedLower;
127 std::unique_ptr<CuMat> m_gpuMatrixReorderedUpper;
129 std::unique_ptr<GpuVector<field_type>> m_gpuMatrixReorderedDiag;
131 GpuVector<int> m_gpuNaturalToReorder;
133 GpuVector<int> m_gpuReorderToNatural;
135 GpuVector<field_type> m_gpuDInv;
137 bool m_splitMatrix;
139 bool m_tuneThreadBlockSizes;
142 int m_upperSolveThreadBlockSize = -1;
143 int m_lowerSolveThreadBlockSize = -1;
144 int m_moveThreadBlockSize = -1;
145 int m_DILUFactorizationThreadBlockSize = -1;
146};
147} // end namespace Opm::gpuistl
148
149#endif
Interface class adding the update() method to the preconditioner interface.
Definition: PreconditionerWithUpdate.hpp:32
DILU preconditioner on the GPU.
Definition: GpuDILU.hpp:44
void computeDiagAndMoveReorderedData(int moveThreadBlockSize, int factorizationThreadBlockSize)
Compute the diagonal of the DILU, and update the data of the reordered matrix.
void apply(X &v, const Y &d) override
Apply the preconditoner.
Y range_type
The range type of the preconditioner.
Definition: GpuDILU.hpp:51
virtual bool hasPerfectUpdate() const override
Definition: GpuDILU.hpp:101
GpuSparseMatrix< field_type > CuMat
The GPU matrix type.
Definition: GpuDILU.hpp:55
typename std::remove_const< M >::type matrix_type
The matrix type the preconditioner is for.
Definition: GpuDILU.hpp:47
void pre(X &x, Y &b) override
Prepare the preconditioner.
void post(X &x) override
Post processing.
typename X::field_type field_type
The field type of the preconditioner.
Definition: GpuDILU.hpp:53
static constexpr bool shouldCallPre()
Definition: GpuDILU.hpp:90
Dune::SolverCategory::Category category() const override
Category of the preconditioner (see SolverCategory::Category)
void tuneThreadBlockSizes()
function that will experimentally tune the thread block sizes of the important cuda kernels
X domain_type
The domain type of the preconditioner.
Definition: GpuDILU.hpp:49
GpuDILU(const M &A, bool splitMatrix, bool tuneKernels)
Constructor.
static constexpr bool shouldCallPost()
Definition: GpuDILU.hpp:96
void update() final
Updates the matrix data.
Definition: autotuner.hpp:29