gpu_preconditioner_utils.hpp
Go to the documentation of this file.
1/*
2 Copyright 2025 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_GPUISTL_DETAIL_GPU_PRECONDITIONER_UTILS_HEADER
21#define OPM_GPUISTL_DETAIL_GPU_PRECONDITIONER_UTILS_HEADER
22
24
25#include <dune/istl/bcrsmatrix.hh>
26
27
28namespace Opm::gpuistl::detail {
29
46template<class Operator, class BlockType>
47Dune::BCRSMatrix<BlockType> makeCPUMatrix(const Operator& op) {
48 // TODO: Make this more efficient. Maybe we can simply copy the memory areas directly?
49 // Do note that this function is anyway going away when we have a GPU
50 // constructor for the preconditioners, so it is not a priority.
51 const auto& gpuMatrix = op.getmat();
52
53 const auto nonZeros = gpuMatrix.getNonZeroValues().asStdVector();
54 const auto rowIndices = gpuMatrix.getRowIndices().asStdVector();
55 const auto columnIndices = gpuMatrix.getColumnIndices().asStdVector();
56
57 const auto numberOfNonZeroes = gpuMatrix.nonzeroes();
58 const auto N = gpuMatrix.N();
59
60 Dune::BCRSMatrix<BlockType> matrix(N, N, numberOfNonZeroes, Dune::BCRSMatrix<BlockType>::row_wise);
61 for (auto row = matrix.createbegin(); row != matrix.createend(); ++row) {
62 for (auto j = rowIndices[row.index()]; j != rowIndices[row.index() + 1]; ++j) {
63 const auto columnIndex = columnIndices[j];
64 row.insert(columnIndex);
65 }
66 }
67
68 for (std::size_t i = 0; i < N; ++i) {
69 for (auto j = rowIndices[i]; j != rowIndices[i + 1]; ++j) {
70 const auto columnIndex = columnIndices[j];
71 // Now it gets a bit tricky, first we need to fetch the block matrix
72 BlockType blockMatrix;
73 constexpr static auto rows = BlockType::rows;
74
75 for (std::size_t k = 0; k < rows; ++k) {
76 for (std::size_t l = 0; l < rows; ++l) {
77 blockMatrix[k][l] = nonZeros[j * rows * rows + k * rows + l];
78 }
79 }
80 matrix[i][columnIndex] = blockMatrix;
81 }
82 }
83
84 return matrix;
85}
86
87
88} // namespace Opm::gpuistl::detail
89
90#endif // OPM_GPUISTL_DETAIL_GPU_PRECONDITIONER_UTILS_HEADER
Definition: autotuner.hpp:30
Dune::BCRSMatrix< BlockType > makeCPUMatrix(const Operator &op)
Utility functions for GPU preconditioner creation.
Definition: gpu_preconditioner_utils.hpp:47