Arcane  v4.1.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
relaxation_cusparse_ilu0.h
1// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
2//-----------------------------------------------------------------------------
3// Copyright 2000-2026 CEA (www.cea.fr) IFPEN (www.ifpenergiesnouvelles.com)
4// See the top-level COPYRIGHT file for details.
5// SPDX-License-Identifier: Apache-2.0
6//-----------------------------------------------------------------------------
7/*---------------------------------------------------------------------------*/
8/* relaxation_cusparse_ilu0.h (C) 2000-2026 */
9/* */
10/* Implementation of ILU0 smoother for CUDA backend. */
11/*---------------------------------------------------------------------------*/
12#ifndef ARCCORE_ALINA_RELAXATION_CUSPARSEILU0_H
13#define ARCCORE_ALINA_RELAXATION_CUSPARSEILU0_H
14/*---------------------------------------------------------------------------*/
15/*---------------------------------------------------------------------------*/
16/*
17 * This file is based on the work on AMGCL library (version march 2026)
18 * which can be found at https://github.com/ddemidov/amgcl.
19 *
20 * Copyright (c) 2012-2022 Denis Demidov <dennis.demidov@gmail.com>
21 * SPDX-License-Identifier: MIT
22 */
23/*---------------------------------------------------------------------------*/
24/*---------------------------------------------------------------------------*/
25
26#include <type_traits>
27
28#include <thrust/device_vector.h>
29#include <cusparse_v2.h>
30
31#include "arccore/alina/CudaBackend.h"
32
33/*---------------------------------------------------------------------------*/
34/*---------------------------------------------------------------------------*/
35
36namespace Arcane::Alina::relaxation
37{
38
39/*---------------------------------------------------------------------------*/
40/*---------------------------------------------------------------------------*/
41
42template <class Backend> struct ilu0;
43
44/*---------------------------------------------------------------------------*/
45/*---------------------------------------------------------------------------*/
49template <typename real>
50struct ilu0<backend::cuda<real>>
51{
52 typedef real value_type;
53 typedef backend::cuda<real> Backend;
54
55 struct params
56 {
58 float damping = 1.0;
59
60 params() = default;
61
62 params(const PropertyTree& p)
63 : ARCCORE_ALINA_PARAMS_IMPORT_VALUE(p, damping)
64 {
65 p.check_params({ "damping" });
66 }
67
68 void get(Alina::PropertyTree& p, const std::string& path) const
69 {
70 ARCCORE_ALINA_PARAMS_EXPORT_VALUE(p, path, damping);
71 }
72 } prm;
73
74 template <class Matrix>
75 ilu0(const Matrix& A, const params& prm, const typename Backend::params& bprm)
76 : prm(prm)
77 , handle(bprm.cusparse_handle)
78 , n(backend::nbRow(A))
79 , nnz(backend::nonzeros(A))
80 , ptr(A.ptr, A.ptr + n + 1)
81 , col(A.col, A.col + nnz)
82 , val(A.val, A.val + nnz)
83 , y(n)
84 {
85 // LU decomposition
86 std::shared_ptr<std::remove_pointer<cusparseMatDescr_t>::type> descr_M;
87 std::shared_ptr<std::remove_pointer<csrilu02Info_t>::type> info_M;
88
89 {
90 cusparseMatDescr_t descr;
91 csrilu02Info_t info;
92
93 ARCCORE_ALINA_CALL_CUDA(cusparseCreateMatDescr(&descr));
94 ARCCORE_ALINA_CALL_CUDA(cusparseSetMatIndexBase(descr, CUSPARSE_INDEX_BASE_ZERO));
95 ARCCORE_ALINA_CALL_CUDA(cusparseSetMatType(descr, CUSPARSE_MATRIX_TYPE_GENERAL));
96
97 ARCCORE_ALINA_CALL_CUDA(cusparseCreateCsrilu02Info(&info));
98
99 descr_M.reset(descr, backend::detail::cuda_deleter());
100 info_M.reset(info, backend::detail::cuda_deleter());
101
102 int buf_size;
103
104 ARCCORE_ALINA_CALL_CUDA(
105 cusparseXcsrilu02_bufferSize(handle, n, nnz, descr_M.get(),
106 thrust::raw_pointer_cast(&val[0]),
107 thrust::raw_pointer_cast(&ptr[0]),
108 thrust::raw_pointer_cast(&col[0]),
109 info_M.get(), &buf_size));
110
111 thrust::device_vector<char> bufLU(buf_size);
112
113 // Analysis and incomplete factorization of the system matrix.
114 int structural_zero;
115 int numerical_zero;
116
117 ARCCORE_ALINA_CALL_CUDA(
118 cusparseXcsrilu02_analysis(handle,
119 n,
120 nnz,
121 descr_M.get(),
122 thrust::raw_pointer_cast(&val[0]),
123 thrust::raw_pointer_cast(&ptr[0]),
124 thrust::raw_pointer_cast(&col[0]),
125 info_M.get(),
126 CUSPARSE_SOLVE_POLICY_USE_LEVEL,
127 thrust::raw_pointer_cast(&bufLU[0])));
128
129 precondition(
130 CUSPARSE_STATUS_ZERO_PIVOT != cusparseXcsrilu02_zeroPivot(handle, info_M.get(), &structural_zero),
131 "Zero pivot in cuSPARSE ILU0");
132
133 ARCCORE_ALINA_CALL_CUDA(
134 cusparseXcsrilu02(handle,
135 n,
136 nnz,
137 descr_M.get(),
138 thrust::raw_pointer_cast(&val[0]),
139 thrust::raw_pointer_cast(&ptr[0]),
140 thrust::raw_pointer_cast(&col[0]),
141 info_M.get(),
142 CUSPARSE_SOLVE_POLICY_USE_LEVEL,
143 thrust::raw_pointer_cast(&bufLU[0])));
144 precondition(
145 CUSPARSE_STATUS_ZERO_PIVOT != cusparseXcsrilu02_zeroPivot(handle, info_M.get(), &numerical_zero),
146 "Zero pivot in cuSPARSE ILU0");
147 }
148
149 // Triangular solvers
150#if CUDART_VERSION >= 11000
151 const real alpha = 1;
152 thrust::device_vector<value_type> t(n);
153
154 descr_y.reset(
155 backend::detail::cuda_vector_description(y),
156 backend::detail::cuda_deleter());
157
158 std::shared_ptr<std::remove_pointer<cusparseDnVecDescr_t>::type> descr_t(
159 backend::detail::cuda_vector_description(t),
160 backend::detail::cuda_deleter());
161
162 cusparseFillMode_t fill_lower = CUSPARSE_FILL_MODE_LOWER;
163 cusparseFillMode_t fill_upper = CUSPARSE_FILL_MODE_UPPER;
164 cusparseDiagType_t diag_unit = CUSPARSE_DIAG_TYPE_UNIT;
165 cusparseDiagType_t diag_non_unit = CUSPARSE_DIAG_TYPE_NON_UNIT;
166
167 // Triangular solver for L
168 {
169 descr_L.reset(
170 backend::detail::cuda_matrix_description(n, n, nnz, ptr, col, val),
171 backend::detail::cuda_deleter());
172
173 ARCCORE_ALINA_CALL_CUDA(
174 cusparseSpMatSetAttribute(descr_L.get(),
175 CUSPARSE_SPMAT_FILL_MODE,
176 &fill_lower,
177 sizeof(fill_lower)));
178
179 ARCCORE_ALINA_CALL_CUDA(
180 cusparseSpMatSetAttribute(descr_L.get(),
181 CUSPARSE_SPMAT_DIAG_TYPE,
182 &diag_unit,
183 sizeof(diag_unit)));
184
185 size_t buf_size;
186
187 cusparseSpSVDescr_t desc;
188 ARCCORE_ALINA_CALL_CUDA(cusparseSpSV_createDescr(&desc));
189 descr_SL.reset(desc, backend::detail::cuda_deleter());
190
191 ARCCORE_ALINA_CALL_CUDA(
192 cusparseSpSV_bufferSize(handle,
193 CUSPARSE_OPERATION_NON_TRANSPOSE,
194 &alpha,
195 descr_L.get(),
196 descr_t.get(),
197 descr_y.get(),
198 backend::detail::cuda_datatype<real>(),
199 CUSPARSE_SPSV_ALG_DEFAULT,
200 descr_SL.get(),
201 &buf_size));
202
203 bufL.resize(buf_size);
204
205 ARCCORE_ALINA_CALL_CUDA(
206 cusparseSpSV_analysis(handle,
207 CUSPARSE_OPERATION_NON_TRANSPOSE,
208 &alpha,
209 descr_L.get(),
210 descr_t.get(),
211 descr_y.get(),
212 backend::detail::cuda_datatype<real>(),
213 CUSPARSE_SPSV_ALG_DEFAULT,
214 descr_SL.get(),
215 thrust::raw_pointer_cast(&bufL[0])));
216 }
217
218 // Triangular solver for U
219 {
220 descr_U.reset(
221 backend::detail::cuda_matrix_description(n, n, nnz, ptr, col, val),
222 backend::detail::cuda_deleter());
223
224 ARCCORE_ALINA_CALL_CUDA(
225 cusparseSpMatSetAttribute(descr_U.get(),
226 CUSPARSE_SPMAT_FILL_MODE,
227 &fill_upper,
228 sizeof(fill_upper)));
229
230 ARCCORE_ALINA_CALL_CUDA(
231 cusparseSpMatSetAttribute(descr_U.get(),
232 CUSPARSE_SPMAT_DIAG_TYPE,
233 &diag_non_unit,
234 sizeof(diag_non_unit)));
235
236 size_t buf_size;
237
238 cusparseSpSVDescr_t desc;
239 ARCCORE_ALINA_CALL_CUDA(cusparseSpSV_createDescr(&desc));
240 descr_SU.reset(desc, backend::detail::cuda_deleter());
241
242 ARCCORE_ALINA_CALL_CUDA(
243 cusparseSpSV_bufferSize(handle,
244 CUSPARSE_OPERATION_NON_TRANSPOSE,
245 &alpha,
246 descr_U.get(),
247 descr_y.get(),
248 descr_t.get(),
249 backend::detail::cuda_datatype<real>(),
250 CUSPARSE_SPSV_ALG_DEFAULT,
251 descr_SU.get(),
252 &buf_size));
253
254 bufU.resize(buf_size);
255
256 ARCCORE_ALINA_CALL_CUDA(
257 cusparseSpSV_analysis(handle,
258 CUSPARSE_OPERATION_NON_TRANSPOSE,
259 &alpha,
260 descr_U.get(),
261 descr_y.get(),
262 descr_t.get(),
263 backend::detail::cuda_datatype<real>(),
264 CUSPARSE_SPSV_ALG_DEFAULT,
265 descr_SU.get(),
266 thrust::raw_pointer_cast(&bufU[0])));
267 }
268#else // CUDART_VERSION >= 11000
269 {
270 cusparseMatDescr_t descr;
271
272 ARCCORE_ALINA_CALL_CUDA(cusparseCreateMatDescr(&descr));
273 ARCCORE_ALINA_CALL_CUDA(cusparseSetMatIndexBase(descr, CUSPARSE_INDEX_BASE_ZERO));
274 ARCCORE_ALINA_CALL_CUDA(cusparseSetMatType(descr, CUSPARSE_MATRIX_TYPE_GENERAL));
275 ARCCORE_ALINA_CALL_CUDA(cusparseSetMatFillMode(descr, CUSPARSE_FILL_MODE_LOWER));
276 ARCCORE_ALINA_CALL_CUDA(cusparseSetMatDiagType(descr, CUSPARSE_DIAG_TYPE_UNIT));
277
278 descr_L.reset(descr, backend::detail::cuda_deleter());
279 }
280 {
281 cusparseMatDescr_t descr;
282
283 ARCCORE_ALINA_CALL_CUDA(cusparseCreateMatDescr(&descr));
284 ARCCORE_ALINA_CALL_CUDA(cusparseSetMatIndexBase(descr, CUSPARSE_INDEX_BASE_ZERO));
285 ARCCORE_ALINA_CALL_CUDA(cusparseSetMatType(descr, CUSPARSE_MATRIX_TYPE_GENERAL));
286 ARCCORE_ALINA_CALL_CUDA(cusparseSetMatFillMode(descr, CUSPARSE_FILL_MODE_UPPER));
287 ARCCORE_ALINA_CALL_CUDA(cusparseSetMatDiagType(descr, CUSPARSE_DIAG_TYPE_NON_UNIT));
288
289 descr_U.reset(descr, backend::detail::cuda_deleter());
290 }
291
292 // Create info structures.
293 {
294 csrsv2Info_t info;
295 ARCCORE_ALINA_CALL_CUDA(cusparseCreateCsrsv2Info(&info));
296 info_L.reset(info, backend::detail::cuda_deleter());
297 }
298 {
299 csrsv2Info_t info;
300 ARCCORE_ALINA_CALL_CUDA(cusparseCreateCsrsv2Info(&info));
301 info_U.reset(info, backend::detail::cuda_deleter());
302 }
303
304 // Allocate scratch buffer.
305 {
306 int buf_size_L;
307 int buf_size_U;
308
309 ARCCORE_ALINA_CALL_CUDA(
310 cusparseXcsrsv2_bufferSize(handle,
311 CUSPARSE_OPERATION_NON_TRANSPOSE,
312 n,
313 nnz,
314 descr_L.get(),
315 thrust::raw_pointer_cast(&val[0]),
316 thrust::raw_pointer_cast(&ptr[0]),
317 thrust::raw_pointer_cast(&col[0]),
318 info_L.get(), &buf_size_L));
319
320 ARCCORE_ALINA_CALL_CUDA(
321 cusparseXcsrsv2_bufferSize(handle,
322 CUSPARSE_OPERATION_NON_TRANSPOSE,
323 n,
324 nnz,
325 descr_U.get(),
326 thrust::raw_pointer_cast(&val[0]),
327 thrust::raw_pointer_cast(&ptr[0]),
328 thrust::raw_pointer_cast(&col[0]),
329 info_U.get(), &buf_size_U));
330
331 buf.resize(std::max(buf_size_L, buf_size_U));
332 }
333
334 ARCCORE_ALINA_CALL_CUDA(
335 cusparseXcsrsv2_analysis(handle,
336 CUSPARSE_OPERATION_NON_TRANSPOSE,
337 n,
338 nnz,
339 descr_L.get(),
340 thrust::raw_pointer_cast(&val[0]),
341 thrust::raw_pointer_cast(&ptr[0]),
342 thrust::raw_pointer_cast(&col[0]),
343 info_L.get(), CUSPARSE_SOLVE_POLICY_USE_LEVEL,
344 thrust::raw_pointer_cast(&buf[0])));
345
346 ARCCORE_ALINA_CALL_CUDA(
347 cusparseXcsrsv2_analysis(handle,
348 CUSPARSE_OPERATION_NON_TRANSPOSE,
349 n,
350 nnz,
351 descr_U.get(),
352 thrust::raw_pointer_cast(&val[0]),
353 thrust::raw_pointer_cast(&ptr[0]),
354 thrust::raw_pointer_cast(&col[0]),
355 info_U.get(), CUSPARSE_SOLVE_POLICY_USE_LEVEL,
356 thrust::raw_pointer_cast(&buf[0])));
357#endif // CUDART_VERSION >= 11000
358 }
359
360 template <class Matrix, class VectorRHS, class VectorX, class VectorTMP>
361 void apply_pre(const Matrix& A, const VectorRHS& rhs, VectorX& x, VectorTMP& tmp) const
362 {
363 backend::residual(rhs, A, x, tmp);
364 solve(tmp);
365 backend::axpby(prm.damping, tmp, 1, x);
366 }
367
368 template <class Matrix, class VectorRHS, class VectorX, class VectorTMP>
369 void apply_post(const Matrix& A, const VectorRHS& rhs, VectorX& x, VectorTMP& tmp) const
370 {
371 backend::residual(rhs, A, x, tmp);
372 solve(tmp);
373 backend::axpby(prm.damping, tmp, 1, x);
374 }
375
376 template <class Matrix, class VectorRHS, class VectorX>
377 void apply(const Matrix& A, const VectorRHS& rhs, VectorX& x) const
378 {
379 backend::copy(rhs, x);
380 solve(x);
381 }
382
383 size_t bytes() const
384 {
385 // This is incomplete, as cusparse structs are opaque.
386 return backend::bytes(ptr) +
387 backend::bytes(col) +
388 backend::bytes(val) +
389 backend::bytes(y) +
390#if CUDART_VERSION >= 11000
391 backend::bytes(bufL) +
392 backend::bytes(bufU)
393#else
394 backend::bytes(buf)
395#endif
396 ;
397 }
398
399 private:
400
401 cusparseHandle_t handle;
402 int n, nnz;
403
404 thrust::device_vector<int> ptr, col;
405 thrust::device_vector<value_type> val;
406 mutable thrust::device_vector<value_type> y;
407
408#if CUDART_VERSION >= 11000
409 std::shared_ptr<std::remove_pointer<cusparseSpMatDescr_t>::type> descr_L, descr_U;
410 std::shared_ptr<std::remove_pointer<cusparseSpSVDescr_t>::type> descr_SL, descr_SU;
411 std::shared_ptr<std::remove_pointer<cusparseDnVecDescr_t>::type> descr_y;
412 mutable thrust::device_vector<char> bufL, bufU;
413#else
414 std::shared_ptr<std::remove_pointer<cusparseMatDescr_t>::type> descr_L, descr_U;
415 std::shared_ptr<std::remove_pointer<csrsv2Info_t>::type> info_L, info_U;
416 mutable thrust::device_vector<char> buf;
417#endif
418
419 template <class VectorX>
420 void solve(VectorX& x) const
421 {
422 value_type alpha = 1;
423
424#if CUDART_VERSION >= 11000
425 std::shared_ptr<std::remove_pointer<cusparseDnVecDescr_t>::type> descr_x(
426 backend::detail::cuda_vector_description(x),
427 backend::detail::cuda_deleter());
428
429 // Solve L * y = x
430 ARCCORE_ALINA_CALL_CUDA(
431 cusparseSpSV_solve(handle,
432 CUSPARSE_OPERATION_NON_TRANSPOSE,
433 &alpha,
434 descr_L.get(),
435 descr_x.get(),
436 descr_y.get(),
437 backend::detail::cuda_datatype<real>(),
438 CUSPARSE_SPSV_ALG_DEFAULT,
439 descr_SL.get()));
440
441 // Solve U * x = y
442 ARCCORE_ALINA_CALL_CUDA(
443 cusparseSpSV_solve(handle,
444 CUSPARSE_OPERATION_NON_TRANSPOSE,
445 &alpha,
446 descr_U.get(),
447 descr_y.get(),
448 descr_x.get(),
449 backend::detail::cuda_datatype<real>(),
450 CUSPARSE_SPSV_ALG_DEFAULT,
451 descr_SU.get()));
452#else // CUDART_VERSION >= 11000
453 // Solve L * y = x
454 ARCCORE_ALINA_CALL_CUDA(
455 cusparseXcsrsv2_solve(handle,
456 CUSPARSE_OPERATION_NON_TRANSPOSE,
457 n,
458 nnz,
459 &alpha,
460 descr_L.get(),
461 thrust::raw_pointer_cast(&val[0]),
462 thrust::raw_pointer_cast(&ptr[0]),
463 thrust::raw_pointer_cast(&col[0]),
464 info_L.get(),
465 thrust::raw_pointer_cast(&x[0]),
466 thrust::raw_pointer_cast(&y[0]),
467 CUSPARSE_SOLVE_POLICY_USE_LEVEL,
468 thrust::raw_pointer_cast(&buf[0])));
469
470 // Solve U * x = y
471 ARCCORE_ALINA_CALL_CUDA(
472 cusparseXcsrsv2_solve(handle,
473 CUSPARSE_OPERATION_NON_TRANSPOSE,
474 n,
475 nnz,
476 &alpha,
477 descr_U.get(),
478 thrust::raw_pointer_cast(&val[0]),
479 thrust::raw_pointer_cast(&ptr[0]),
480 thrust::raw_pointer_cast(&col[0]),
481 info_U.get(),
482 thrust::raw_pointer_cast(&y[0]),
483 thrust::raw_pointer_cast(&x[0]),
484 CUSPARSE_SOLVE_POLICY_USE_LEVEL,
485 thrust::raw_pointer_cast(&buf[0])));
486#endif // CUDART_VERSION >= 11000
487 }
488
489 static cusparseStatus_t
490 cusparseXcsrilu02_bufferSize(cusparseHandle_t handle,
491 int m,
492 int nnz,
493 const cusparseMatDescr_t descrA,
494 double* csrSortedValA,
495 const int* csrSortedRowPtrA,
496 const int* csrSortedColIndA,
497 csrilu02Info_t info,
498 int* pBufferSizeInBytes)
499 {
500 return cusparseDcsrilu02_bufferSize(handle, m, nnz, descrA, csrSortedValA, csrSortedRowPtrA,
501 csrSortedColIndA, info, pBufferSizeInBytes);
502 }
503
504 static cusparseStatus_t
505 cusparseXcsrilu02_bufferSize(cusparseHandle_t handle,
506 int m,
507 int nnz,
508 const cusparseMatDescr_t descrA,
509 float* csrSortedValA,
510 const int* csrSortedRowPtrA,
511 const int* csrSortedColIndA,
512 csrilu02Info_t info,
513 int* pBufferSizeInBytes)
514 {
515 return cusparseScsrilu02_bufferSize(handle, m, nnz, descrA, csrSortedValA, csrSortedRowPtrA,
516 csrSortedColIndA, info, pBufferSizeInBytes);
517 }
518
519 static cusparseStatus_t
520 cusparseXcsrilu02_analysis(cusparseHandle_t handle,
521 int m,
522 int nnz,
523 const cusparseMatDescr_t descrA,
524 const double* csrSortedValA,
525 const int* csrSortedRowPtrA,
526 const int* csrSortedColIndA,
527 csrilu02Info_t info,
528 cusparseSolvePolicy_t policy,
529 void* pBuffer)
530 {
531 return cusparseDcsrilu02_analysis(handle, m, nnz, descrA, csrSortedValA,
532 csrSortedRowPtrA, csrSortedColIndA, info, policy, pBuffer);
533 }
534
535 static cusparseStatus_t
536 cusparseXcsrilu02_analysis(cusparseHandle_t handle,
537 int m,
538 int nnz,
539 const cusparseMatDescr_t descrA,
540 const float* csrSortedValA,
541 const int* csrSortedRowPtrA,
542 const int* csrSortedColIndA,
543 csrilu02Info_t info,
544 cusparseSolvePolicy_t policy,
545 void* pBuffer)
546 {
547 return cusparseScsrilu02_analysis(handle, m, nnz, descrA, csrSortedValA,
548 csrSortedRowPtrA, csrSortedColIndA, info, policy, pBuffer);
549 }
550
551 static cusparseStatus_t
552 cusparseXcsrilu02(cusparseHandle_t handle,
553 int m,
554 int nnz,
555 const cusparseMatDescr_t descrA,
556 double* csrSortedValA_valM,
557 const int* csrSortedRowPtrA,
558 const int* csrSortedColIndA,
559 csrilu02Info_t info,
560 cusparseSolvePolicy_t policy,
561 void* pBuffer)
562 {
563 return cusparseDcsrilu02(handle, m, nnz, descrA,
564 csrSortedValA_valM, csrSortedRowPtrA, csrSortedColIndA,
565 info, policy, pBuffer);
566 }
567
568 static cusparseStatus_t
569 cusparseXcsrilu02(cusparseHandle_t handle,
570 int m,
571 int nnz,
572 const cusparseMatDescr_t descrA,
573 float* csrSortedValA_valM,
574 const int* csrSortedRowPtrA,
575 const int* csrSortedColIndA,
576 csrilu02Info_t info,
577 cusparseSolvePolicy_t policy,
578 void* pBuffer)
579 {
580 return cusparseScsrilu02(handle, m, nnz, descrA,
581 csrSortedValA_valM, csrSortedRowPtrA, csrSortedColIndA,
582 info, policy, pBuffer);
583 }
584
585#if CUDART_VERSION < 11000
586 static cusparseStatus_t
587 cusparseXcsrsv2_bufferSize(cusparseHandle_t handle,
588 cusparseOperation_t transA,
589 int m,
590 int nnz,
591 const cusparseMatDescr_t descrA,
592 double* csrSortedValA,
593 const int* csrSortedRowPtrA,
594 const int* csrSortedColIndA,
595 csrsv2Info_t info,
596 int* pBufferSizeInBytes)
597 {
598 return cusparseDcsrsv2_bufferSize(handle, transA, m, nnz, descrA, csrSortedValA,
599 csrSortedRowPtrA, csrSortedColIndA, info, pBufferSizeInBytes);
600 }
601
602 static cusparseStatus_t
603 cusparseXcsrsv2_bufferSize(cusparseHandle_t handle,
604 cusparseOperation_t transA,
605 int m,
606 int nnz,
607 const cusparseMatDescr_t descrA,
608 float* csrSortedValA,
609 const int* csrSortedRowPtrA,
610 const int* csrSortedColIndA,
611 csrsv2Info_t info,
612 int* pBufferSizeInBytes)
613 {
614 return cusparseScsrsv2_bufferSize(handle, transA, m, nnz, descrA, csrSortedValA,
615 csrSortedRowPtrA, csrSortedColIndA, info, pBufferSizeInBytes);
616 }
617
618 static cusparseStatus_t
619 cusparseXcsrsv2_analysis(cusparseHandle_t handle,
620 cusparseOperation_t transA,
621 int m,
622 int nnz,
623 const cusparseMatDescr_t descrA,
624 const double* csrSortedValA,
625 const int* csrSortedRowPtrA,
626 const int* csrSortedColIndA,
627 csrsv2Info_t info,
628 cusparseSolvePolicy_t policy,
629 void* pBuffer)
630 {
631 return cusparseDcsrsv2_analysis(handle, transA, m, nnz, descrA, csrSortedValA,
632 csrSortedRowPtrA, csrSortedColIndA, info, policy, pBuffer);
633 }
634
635 static cusparseStatus_t
636 cusparseXcsrsv2_analysis(cusparseHandle_t handle,
637 cusparseOperation_t transA,
638 int m,
639 int nnz,
640 const cusparseMatDescr_t descrA,
641 const float* csrSortedValA,
642 const int* csrSortedRowPtrA,
643 const int* csrSortedColIndA,
644 csrsv2Info_t info,
645 cusparseSolvePolicy_t policy,
646 void* pBuffer)
647 {
648 return cusparseScsrsv2_analysis(handle, transA, m, nnz, descrA, csrSortedValA,
649 csrSortedRowPtrA, csrSortedColIndA, info, policy, pBuffer);
650 }
651
652 static cusparseStatus_t
653 cusparseXcsrsv2_solve(cusparseHandle_t handle,
654 cusparseOperation_t transA,
655 int m,
656 int nnz,
657 const double* alpha,
658 const cusparseMatDescr_t descrA,
659 const double* csrSortedValA,
660 const int* csrSortedRowPtrA,
661 const int* csrSortedColIndA,
662 csrsv2Info_t info,
663 const double* f,
664 double* x,
665 cusparseSolvePolicy_t policy,
666 void* pBuffer)
667 {
668 return cusparseDcsrsv2_solve(handle, transA, m,
669 nnz, alpha, descrA, csrSortedValA, csrSortedRowPtrA,
670 csrSortedColIndA, info, f, x, policy, pBuffer);
671 }
672
673 static cusparseStatus_t
674 cusparseXcsrsv2_solve(cusparseHandle_t handle,
675 cusparseOperation_t transA,
676 int m,
677 int nnz,
678 const float* alpha,
679 const cusparseMatDescr_t descrA,
680 const float* csrSortedValA,
681 const int* csrSortedRowPtrA,
682 const int* csrSortedColIndA,
683 csrsv2Info_t info,
684 const float* f,
685 float* x,
686 cusparseSolvePolicy_t policy,
687 void* pBuffer)
688 {
689 return cusparseScsrsv2_solve(handle, transA, m,
690 nnz, alpha, descrA, csrSortedValA, csrSortedRowPtrA,
691 csrSortedColIndA, info, f, x, policy, pBuffer);
692 }
693#endif // CUDART_VERSION < 11000
694};
695
696/*---------------------------------------------------------------------------*/
697/*---------------------------------------------------------------------------*/
698
699} // namespace Arcane::Alina::relaxation
700
701/*---------------------------------------------------------------------------*/
702/*---------------------------------------------------------------------------*/
703
704#endif
Matrix class, to be used by user.