gpusparse_matrix_utilities.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#ifndef OPM_GPUISTL_GPUSPARSE_MATRIX_UTILITIES_HPP
20#define OPM_GPUISTL_GPUSPARSE_MATRIX_UTILITIES_HPP
21
22#include <opm/common/ErrorMacros.hpp>
25
26#include <fmt/core.h>
27#include <memory>
28#include <vector>
29
30#include <cusparse.h>
31
33{
34
39inline auto
41{
42 auto deleter = [](cusparseDnVecDescr_t* desc) {
43 if (desc && *desc) {
44 OPM_CUSPARSE_WARN_IF_ERROR(cusparseDestroyDnVec(*desc));
45 }
46 delete desc;
47 };
48 return std::unique_ptr<cusparseDnVecDescr_t, decltype(deleter)>(new cusparseDnVecDescr_t {nullptr}, deleter);
49}
50
55inline auto
57{
58 auto deleter = [](cusparseSpMatDescr_t* desc) {
59 if (desc && *desc) {
60 OPM_CUSPARSE_WARN_IF_ERROR(cusparseDestroySpMat(*desc));
61 }
62 delete desc;
63 };
64 return std::unique_ptr<cusparseSpMatDescr_t, decltype(deleter)>(new cusparseSpMatDescr_t {nullptr}, deleter);
65}
66
78template <class T, class M>
79std::vector<T>
80extractNonzeroValues(const M& matrix)
81{
82 const std::size_t blockSize = matrix[0][0].N();
83 const std::size_t numberOfNonzeroBlocks = matrix.nonzeroes();
84 const std::size_t numberOfNonzeroElements = blockSize * blockSize * numberOfNonzeroBlocks;
85
86 std::vector<T> nonZeroElementsData(numberOfNonzeroElements);
87
88 std::size_t ind = 0;
89 for (auto& row : matrix) {
90 for (const auto& element : row) {
91 for (std::size_t c = 0; c < blockSize; ++c) {
92 for (std::size_t d = 0; d < blockSize; ++d) {
93 nonZeroElementsData[ind++] = element[c][d];
94 }
95 }
96 }
97 }
98
99 return nonZeroElementsData;
100}
101
112template <class M>
113std::pair<std::vector<int>, std::vector<int>>
115{
116 std::vector<int> columnIndices(matrix.nonzeroes());
117 std::vector<int> rowIndices(matrix.N() + 1);
118
119 rowIndices[0] = 0;
120 std::size_t colIndex = 0;
121 std::size_t rowIndex = 1;
122 for (const auto& row : matrix) {
123 for (auto columnIterator = row.begin(); columnIterator != row.end(); ++columnIterator) {
124 columnIndices[colIndex++] = columnIterator.index();
125 }
126 rowIndices[rowIndex++] = colIndex;
127 }
128
129 return {columnIndices, rowIndices};
130}
131
143template <class T, class M>
144T*
146{
147 // We must find the pointer to the first element in the matrix
148 // Iterate until we find an element, we can get the blocksize from the element
149 // TODO: Can this be done more cleanly in the DUNE api to access the raw data more directly?
150 T* nonZeroElementsTmp = nullptr;
151 for (auto rowIt = matrix.begin(); rowIt != matrix.end(); ++rowIt) {
152 auto colIt = rowIt->begin();
153 if (colIt != rowIt->end()) {
154 nonZeroElementsTmp = const_cast<T*>(&((*colIt)[0][0]));
155 break;
156 }
157 }
158 OPM_ERROR_IF(nonZeroElementsTmp == nullptr, "No non-zero elements found in matrix");
159 return nonZeroElementsTmp;
160}
161
173template <class M>
174void
175validateMatrixConversion(const M& matrix, const std::vector<int>& columnIndices, const std::vector<int>& rowIndices)
176{
177 const size_t numberOfRows = matrix.N();
178 const size_t numberOfNonzeroBlocks = matrix.nonzeroes();
179
180 OPM_ERROR_IF(rowIndices[matrix.N()] != to_int(matrix.nonzeroes()),
181 "Error: Row indices do not sum to number of nonzeroes");
182 OPM_ERROR_IF(rowIndices.size() != numberOfRows + 1, "Error: Row indices size does not match matrix dimensions");
183 OPM_ERROR_IF(columnIndices.size() != numberOfNonzeroBlocks,
184 "Error: Column indices size does not match number of nonzero blocks");
185}
186
197inline void
199 size_t nonzeroes1, size_t blockSize1, size_t N1, size_t nonzeroes2, size_t blockSize2, size_t N2)
200{
201 OPM_ERROR_IF(nonzeroes1 != nonzeroes2, "Matrix does not have the same number of non-zero elements.");
202 OPM_ERROR_IF(blockSize1 != blockSize2, "Matrix does not have the same blocksize.");
203 OPM_ERROR_IF(N1 != N2, "Matrix does not have the same number of rows.");
204}
205
213inline void
214validateVectorMatrixSizes(size_t vectorSize, size_t matrixBlockSize, size_t matrixColumns)
215{
216 const size_t expectedSize = matrixBlockSize * matrixColumns;
217 if (vectorSize != expectedSize) {
218 OPM_THROW(std::invalid_argument,
219 fmt::format("Size mismatch. Input vector has {} elements, while matrix expects {} elements.",
220 vectorSize,
221 expectedSize));
222 }
223}
224
231#define INSTANTIATE_SPARSE_MATRIX_DUNE_OPERATIONS(CLASS_NAME, realtype, blockdim) \
232 template CLASS_NAME<realtype> CLASS_NAME<realtype>::fromMatrix( \
233 const Dune::BCRSMatrix<Dune::FieldMatrix<realtype, blockdim, blockdim>>&, bool); \
234 template CLASS_NAME<realtype> CLASS_NAME<realtype>::fromMatrix( \
235 const Dune::BCRSMatrix<Opm::MatrixBlock<realtype, blockdim, blockdim>>&, bool); \
236 template void CLASS_NAME<realtype>::updateNonzeroValues( \
237 const Dune::BCRSMatrix<Dune::FieldMatrix<realtype, blockdim, blockdim>>&, bool); \
238 template void CLASS_NAME<realtype>::updateNonzeroValues( \
239 const Dune::BCRSMatrix<Opm::MatrixBlock<realtype, blockdim, blockdim>>&, bool)
240
244#define INSTANTIATE_FOR_TYPE_AND_CLASS(CLASS_NAME, T) \
245 INSTANTIATE_SPARSE_MATRIX_DUNE_OPERATIONS(CLASS_NAME, T, 1); \
246 INSTANTIATE_SPARSE_MATRIX_DUNE_OPERATIONS(CLASS_NAME, T, 2); \
247 INSTANTIATE_SPARSE_MATRIX_DUNE_OPERATIONS(CLASS_NAME, T, 3); \
248 INSTANTIATE_SPARSE_MATRIX_DUNE_OPERATIONS(CLASS_NAME, T, 4); \
249 INSTANTIATE_SPARSE_MATRIX_DUNE_OPERATIONS(CLASS_NAME, T, 5); \
250 INSTANTIATE_SPARSE_MATRIX_DUNE_OPERATIONS(CLASS_NAME, T, 6)
251
252} // namespace Opm::gpuistl::detail
253
254#endif
#define OPM_CUSPARSE_WARN_IF_ERROR(expression)
OPM_CUSPARSE_WARN_IF_ERROR checks the return type of the cusparse expression (function call) and issu...
Definition: cusparse_safe_call.hpp:206
Definition: autotuner.hpp:30
int to_int(std::size_t s)
to_int converts a (on most relevant platforms) 64 bits unsigned size_t to a signed 32 bits signed int
Definition: safe_conversion.hpp:52
std::vector< T > extractNonzeroValues(const M &matrix)
Extract non-zero values from a DUNE matrix in the order expected by cuSPARSE.
Definition: gpusparse_matrix_utilities.hpp:80
T * findFirstElementPointer(const M &matrix)
Find pointer to first matrix element for direct memory access.
Definition: gpusparse_matrix_utilities.hpp:145
void validateMatrixConversion(const M &matrix, const std::vector< int > &columnIndices, const std::vector< int > &rowIndices)
Validate matrix conversion parameters.
Definition: gpusparse_matrix_utilities.hpp:175
void validateMatrixCompatibility(size_t nonzeroes1, size_t blockSize1, size_t N1, size_t nonzeroes2, size_t blockSize2, size_t N2)
Validate that two matrices have compatible dimensions.
Definition: gpusparse_matrix_utilities.hpp:198
void validateVectorMatrixSizes(size_t vectorSize, size_t matrixBlockSize, size_t matrixColumns)
Validate vector-matrix size compatibility for operations.
Definition: gpusparse_matrix_utilities.hpp:214
auto makeSafeVectorDescriptor()
Create RAII-managed cuSPARSE dense vector descriptor.
Definition: gpusparse_matrix_utilities.hpp:40
std::pair< std::vector< int >, std::vector< int > > extractSparsityPattern(const M &matrix)
Extract sparsity pattern from a DUNE matrix.
Definition: gpusparse_matrix_utilities.hpp:114
auto makeSafeMatrixDescriptor()
Create RAII-managed cuSPARSE sparse matrix descriptor.
Definition: gpusparse_matrix_utilities.hpp:56