cusparse_matrix_operations.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_CUISTL_CUSPARSE_MATRIX_OPERATIONS_HPP
20#define OPM_CUISTL_CUSPARSE_MATRIX_OPERATIONS_HPP
21#include <cstddef>
22#include <vector>
23namespace Opm::cuistl::detail
24{
25
34template <class T, int blocksize>
35void invertDiagonalAndFlatten(T* mat, int* rowIndices, int* colIndices, size_t numberOfRows, T* vec);
36
53template <class T, int blocksize>
54void computeLowerSolveLevelSet(T* reorderedMat,
55 int* rowIndices,
56 int* colIndices,
57 int* indexConversion,
58 int startIdx,
59 int rowsInLevelSet,
60 const T* dInv,
61 const T* d,
62 T* v);
63
79template <class T, int blocksize>
80void computeUpperSolveLevelSet(T* reorderedMat,
81 int* rowIndices,
82 int* colIndices,
83 int* indexConversion,
84 int startIdx,
85 int rowsInLevelSet,
86 const T* dInv,
87 T* v);
88
105template <class T, int blocksize>
106void computeDiluDiagonal(T* reorderedMat,
107 int* rowIndices,
108 int* colIndices,
109 int* reorderedToNatural,
110 int* naturalToReordered,
111 int startIdx,
112 int rowsInLevelSet,
113 T* dInv);
114
125template <class T, int blocksize>
127 T* srcMatrix, int* srcRowIndices, T* dstMatrix, int* dstRowIndices, int* naturalToReordered, size_t numberOfRows);
128} // namespace Opm::cuistl::detail
129#endif
Definition: cublas_safe_call.hpp:32
void computeDiluDiagonal(T *reorderedMat, int *rowIndices, int *colIndices, int *reorderedToNatural, int *naturalToReordered, int startIdx, int rowsInLevelSet, T *dInv)
Computes the ILU0 of the diagonal elements of the reordered matrix and stores it in a reordered vecto...
void computeLowerSolveLevelSet(T *reorderedMat, int *rowIndices, int *colIndices, int *indexConversion, int startIdx, int rowsInLevelSet, const T *dInv, const T *d, T *v)
Perform a lower solve on certain rows in a matrix that can safely be computed in parallel.
void computeUpperSolveLevelSet(T *reorderedMat, int *rowIndices, int *colIndices, int *indexConversion, int startIdx, int rowsInLevelSet, const T *dInv, T *v)
Perform an upper solve on certain rows in a matrix that can safely be computed in parallel.
void invertDiagonalAndFlatten(T *mat, int *rowIndices, int *colIndices, size_t numberOfRows, T *vec)
This function receives a matrix, and the inverse of the matrix containing only its diagonal is stored...
void copyMatDataToReordered(T *srcMatrix, int *srcRowIndices, T *dstMatrix, int *dstRowIndices, int *naturalToReordered, size_t numberOfRows)
Reorders the elements of a matrix by copying them from one matrix to another using a permutation list...