DILUKernels.hpp
Go to the documentation of this file.
1/*
2 Copyright 2024 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_DILU_KERNELS_HPP
20#define OPM_DILU_KERNELS_HPP
21
22#include <cstddef>
23#include <cuda.h>
24#include <cuda_runtime.h>
26#include <vector>
27
29{
30
47template <class T, int blocksize>
48void solveLowerLevelSet(T* reorderedMat,
49 int* rowIndices,
50 int* colIndices,
51 int* indexConversion,
52 int startIdx,
53 int rowsInLevelSet,
54 const T* dInv,
55 const T* d,
56 T* v,
57 int threadBlockSize,
58 cudaStream_t stream);
59
73template <class T, int blocksize>
75 const int* rowIndices,
76 const int* colIndices,
77 const size_t* indexConversion,
78 int startIdx,
79 int rowsInLevelSet,
80 const T* dInv,
81 const T* d,
82 T* v,
83 int threadBlockSize,
84 cudaStream_t stream);
85
102template <int blocksize, class LinearSolverScalar, class MatrixScalar, class DiagonalScalar>
103void solveLowerLevelSetSplit(MatrixScalar* reorderedUpperMat,
104 int* rowIndices,
105 int* colIndices,
106 int* indexConversion,
107 int startIdx,
108 int rowsInLevelSet,
109 const DiagonalScalar* dInv,
110 const LinearSolverScalar* d,
111 LinearSolverScalar* v,
112 int threadBlockSize,
113 cudaStream_t stream);
114
130template <class T, int blocksize>
131void solveUpperLevelSet(T* reorderedMat,
132 int* rowIndices,
133 int* colIndices,
134 int* indexConversion,
135 int startIdx,
136 int rowsInLevelSet,
137 const T* dInv,
138 T* v,
139 int threadBlockSize,
140 cudaStream_t stream);
141
155template <class T, int blocksize>
157 const int* rowIndices,
158 const int* colIndices,
159 const size_t* indexConversion,
160 int startIdx,
161 int rowsInLevelSet,
162 const T* dInv,
163 T* v,
164 int threadBlockSize,
165 cudaStream_t stream);
166
182template <int blocksize, class LinearSolverScalar, class MatrixScalar, class DiagonalScalar>
183void solveUpperLevelSetSplit(MatrixScalar* reorderedUpperMat,
184 int* rowIndices,
185 int* colIndices,
186 int* indexConversion,
187 int startIdx,
188 int rowsInLevelSet,
189 const DiagonalScalar* dInv,
190 LinearSolverScalar* v,
191 int threadBlockSize,
192 cudaStream_t stream);
193
210template <class T, int blocksize>
211void computeDiluDiagonal(T* reorderedMat,
212 int* rowIndices,
213 int* colIndices,
214 int* reorderedToNatural,
215 int* naturalToReordered,
216 size_t* diagIndices,
217 int startIdx,
218 int rowsInLevelSet,
219 T* dInv,
220 int threadBlockSize);
221
234template <class T, int blocksize>
236 const int* rowIndices,
237 const int* colIndices,
238 const size_t* indexConversion,
239 const size_t* diagIndices,
240 int startIdx,
241 int rowsInLevelSet,
242 T* dInv,
243 int threadBlockSize);
244
266template <int blocksize, class InputScalar, class OutputScalar, MatrixStorageMPScheme>
267void computeDiluDiagonalSplit(const InputScalar* srcReorderedLowerMat,
268 int* lowerRowIndices,
269 int* lowerColIndices,
270 const InputScalar* srcReorderedUpperMat,
271 int* upperRowIndices,
272 int* upperColIndices,
273 const InputScalar* srcDiagonal,
274 int* reorderedToNatural,
275 int* naturalToReordered,
276 int startIdx,
277 int rowsInLevelSet,
278 InputScalar* dInv,
279 OutputScalar* dstDiagonal,
280 OutputScalar* dstLowerMat,
281 OutputScalar* dstUpperMat,
282 int threadBlockSize);
283
284} // namespace Opm::gpuistl::detail::DILU
285#endif
Definition: DILUKernels.hpp:29
void computeDiluDiagonal(T *reorderedMat, int *rowIndices, int *colIndices, int *reorderedToNatural, int *naturalToReordered, size_t *diagIndices, int startIdx, int rowsInLevelSet, T *dInv, int threadBlockSize)
Computes the ILU0 of the diagonal elements of the reordered matrix and stores it in a reordered vecto...
void solveUpperLevelSetNoReorder(const T *mat, const int *rowIndices, const int *colIndices, const size_t *indexConversion, int startIdx, int rowsInLevelSet, const T *dInv, T *v, int threadBlockSize, cudaStream_t stream)
Perform an upper solve on certain rows in a matrix that can safely be computed in parallel.
void computeDiluDiagonalNoReorder(const T *mat, const int *rowIndices, const int *colIndices, const size_t *indexConversion, const size_t *diagIndices, int startIdx, int rowsInLevelSet, T *dInv, int threadBlockSize)
Computes the DILU of a BCSR matrix and stores it in a vector containing the diagonal blocks.
void solveLowerLevelSetNoReorder(const T *mat, const int *rowIndices, const int *colIndices, const size_t *indexConversion, int startIdx, int rowsInLevelSet, const T *dInv, const T *d, T *v, int threadBlockSize, cudaStream_t stream)
Perform a lower solve on certain rows in a matrix that can safely be computed in parallel.
void solveLowerLevelSetSplit(MatrixScalar *reorderedUpperMat, int *rowIndices, int *colIndices, int *indexConversion, int startIdx, int rowsInLevelSet, const DiagonalScalar *dInv, const LinearSolverScalar *d, LinearSolverScalar *v, int threadBlockSize, cudaStream_t stream)
Perform a lower solve on certain rows in a matrix that can safely be computed in parallel.
void solveLowerLevelSet(T *reorderedMat, int *rowIndices, int *colIndices, int *indexConversion, int startIdx, int rowsInLevelSet, const T *dInv, const T *d, T *v, int threadBlockSize, cudaStream_t stream)
Perform a lower solve on certain rows in a matrix that can safely be computed in parallel.
void solveUpperLevelSet(T *reorderedMat, int *rowIndices, int *colIndices, int *indexConversion, int startIdx, int rowsInLevelSet, const T *dInv, T *v, int threadBlockSize, cudaStream_t stream)
Perform an upper solve on certain rows in a matrix that can safely be computed in parallel.
void computeDiluDiagonalSplit(const InputScalar *srcReorderedLowerMat, int *lowerRowIndices, int *lowerColIndices, const InputScalar *srcReorderedUpperMat, int *upperRowIndices, int *upperColIndices, const InputScalar *srcDiagonal, int *reorderedToNatural, int *naturalToReordered, int startIdx, int rowsInLevelSet, InputScalar *dInv, OutputScalar *dstDiagonal, OutputScalar *dstLowerMat, OutputScalar *dstUpperMat, int threadBlockSize)
Computes the ILU0 of the diagonal elements of the split reordered matrix and stores it in a reordered...
void solveUpperLevelSetSplit(MatrixScalar *reorderedUpperMat, int *rowIndices, int *colIndices, int *indexConversion, int startIdx, int rowsInLevelSet, const DiagonalScalar *dInv, LinearSolverScalar *v, int threadBlockSize, cudaStream_t stream)
Perform an upper solve on certain rows in a matrix that can safely be computed in parallel.