19#ifndef OPM_GPUISTL_GPUSPARSE_MATRIX_UTILITIES_HPP
20#define OPM_GPUISTL_GPUSPARSE_MATRIX_UTILITIES_HPP
22#include <opm/common/ErrorMacros.hpp>
42 auto deleter = [](cusparseDnVecDescr_t* desc) {
48 return std::unique_ptr<cusparseDnVecDescr_t, decltype(deleter)>(
new cusparseDnVecDescr_t {
nullptr}, deleter);
58 auto deleter = [](cusparseSpMatDescr_t* desc) {
64 return std::unique_ptr<cusparseSpMatDescr_t, decltype(deleter)>(
new cusparseSpMatDescr_t {
nullptr}, deleter);
78template <
class T,
class M>
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;
86 std::vector<T> nonZeroElementsData(numberOfNonzeroElements);
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];
99 return nonZeroElementsData;
113std::pair<std::vector<int>, std::vector<int>>
116 std::vector<int> columnIndices(matrix.nonzeroes());
117 std::vector<int> rowIndices(matrix.N() + 1);
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();
126 rowIndices[rowIndex++] = colIndex;
129 return {columnIndices, rowIndices};
143template <
class T,
class M>
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]));
158 OPM_ERROR_IF(nonZeroElementsTmp ==
nullptr,
"No non-zero elements found in matrix");
159 return nonZeroElementsTmp;
177 const size_t numberOfRows = matrix.N();
178 const size_t numberOfNonzeroBlocks = matrix.nonzeroes();
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");
199 size_t nonzeroes1,
size_t blockSize1,
size_t N1,
size_t nonzeroes2,
size_t blockSize2,
size_t N2)
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.");
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.",
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)
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)
#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