aoclsparse_?mv() - 5.2 English - 68552

AOCL API Guide (68552)

Document ID
68552
Release Date
2025-12-29
Version
5.2 English
aoclsparse_status aoclsparse_smv(aoclsparse_operation op, const float *alpha, aoclsparse_matrix A, const aoclsparse_mat_descr descr, const float *x, const float *beta, float *y)#
aoclsparse_status aoclsparse_dmv(aoclsparse_operation op, const double *alpha, aoclsparse_matrix A, const aoclsparse_mat_descr descr, const double *x, const double *beta, double *y)#
aoclsparse_status aoclsparse_cmv(aoclsparse_operation op, const aoclsparse_float_complex *alpha, aoclsparse_matrix A, const aoclsparse_mat_descr descr, const aoclsparse_float_complex *x, const aoclsparse_float_complex *beta, aoclsparse_float_complex *y)#
aoclsparse_status aoclsparse_zmv(aoclsparse_operation op, const aoclsparse_double_complex *alpha, aoclsparse_matrix A, const aoclsparse_mat_descr descr, const aoclsparse_double_complex *x, const aoclsparse_double_complex *beta, aoclsparse_double_complex *y)#

Compute sparse matrix-vector multiplication for real/complex single and double data precisions.

The aoclsparse_?mv perform sparse matrix-vector products of the form

\[ y = \alpha \, op(A) \, x + \beta \, y, \]
where, \(x\) and \(y\) are dense vectors, \(\alpha\) and \(\beta\) are scalars, and \(A\) is a sparse matrix structure. The matrix operation \(op()\) is defined as:
\[\begin{split} op(A) = \left\{ \begin{array}{ll} A, & \text{if } {\bf\mathsf{op}} = \text{aoclsparse}\_\text{operation}\_\text{none} \\ A^T, & \text{if } {\bf\mathsf{op}} = \text{aoclsparse}\_\text{operation}\_\text{transpose} \\ A^H, & \text{if } {\bf\mathsf{op}} = \text{aoclsparse}\_\text{operation}\_\text{conjugate}\_\text{transpose} \end{array} \right. \end{split}\]

 1/* ************************************************************************
 2 * Copyright (c) 2020-2025 Advanced Micro Devices, Inc.
 3 *
 4 * Permission is hereby granted, free of charge, to any person obtaining a copy
 5 * of this software and associated documentation files (the "Software"), to deal
 6 * in the Software without restriction, including without limitation the rights
 7 * to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
 8 * copies of the Software, and to permit persons to whom the Software is
 9 * furnished to do so, subject to the following conditions:
10 *
11 * The above copyright notice and this permission notice shall be included in
12 * all copies or substantial portions of the Software.
13 *
14 * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
15 * IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
16 * FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
17 * AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
18 * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
19 * OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN
20 * THE SOFTWARE.
21 *
22 * ************************************************************************ */
23/*
24    Description: In this example, we create a sparse matrix in CSR format
25                 and perform a matrix-vector multiplication (SpMV) using the
26                 C++ interface of aoclsparse.
27*/
28#include "aoclsparse.hpp"
29
30#include <iostream>
31#include <vector>
32
33int main(void)
34{
35    aoclsparse_operation trans = aoclsparse_operation_none;
36
37    double alpha = 1.0;
38    double beta  = 0.0;
39
40    const aoclsparse_int m   = 5; // Number of rows
41    const aoclsparse_int n   = 5; // Number of columns
42    const aoclsparse_int nnz = 8; // Number of non-zero elements
43
44    // Print aoclsparse version
45    std::cout << aoclsparse_get_version() << std::endl;
46
47    // Create matrix descriptor
48    aoclsparse_mat_descr descr;
49    // aoclsparse_create_mat_descr set aoclsparse_matrix_type to aoclsparse_matrix_type_general
50    // and aoclsparse_index_base to aoclsparse_index_base_zero.
51    aoclsparse_create_mat_descr(&descr);
52
53    aoclsparse_index_base base = aoclsparse_index_base_zero;
54
55    // Initialise matrix
56    //  1  0  0  2  0
57    //  0  3  0  0  0
58    //  0  0  4  0  0
59    //  0  5  0  6  7
60    //  0  0  0  0  8
61    std::vector<aoclsparse_int> csr_row     = {0, 2, 3, 4, 7, 8};
62    std::vector<aoclsparse_int> csr_col_ind = {0, 3, 1, 2, 1, 3, 4, 4};
63    std::vector<double>         csr_val     = {1, 2, 3, 4, 5, 6, 7, 8};
64    aoclsparse_matrix           A;
65    auto                        status = aoclsparse::create_csr(
66        &A, base, m, n, nnz, csr_row.data(), csr_col_ind.data(), csr_val.data());
67
68    if(status != aoclsparse_status_success)
69    {
70        std::cout << "Error creating the matrix, status = " << status << std::endl;
71        return 1;
72    }
73
74    // Initialise vectors
75    std::vector<double> x = {1.0, 2.0, 3.0, 4.0, 5.0};
76    std::vector<double> y(m, 0.0); // Output vector initialized to zero
77
78    //to identify hint id(which routine is to be executed, destroyed later)
79    aoclsparse_set_mv_hint(A, trans, descr, 1);
80
81    // Optimize the matrix, "A"
82    aoclsparse_optimize(A);
83
84    std::cout << "Invoking aoclsparse_dmv..";
85    //Invoke SPMV API (double precision)
86    aoclsparse::mv(trans, &alpha, A, descr, x.data(), &beta, y.data());
87
88    std::cout << "Done." << std::endl;
89
90    std::cout << "Output Vector:" << std::endl;
91    for(aoclsparse_int i = 0; i < m; i++)
92        std::cout << y[i] << std::endl;
93
94    aoclsparse_destroy_mat_descr(descr);
95    aoclsparse_destroy(&A);
96    return 0;
97}

  1/* ************************************************************************
  2 * Copyright (c) 2023 Advanced Micro Devices, Inc.
  3 *
  4 * Permission is hereby granted, free of charge, to any person obtaining a copy
  5 * of this software and associated documentation files (the "Software"), to deal
  6 * in the Software without restriction, including without limitation the rights
  7 * to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
  8 * copies of the Software, and to permit persons to whom the Software is
  9 * furnished to do so, subject to the following conditions:
 10 *
 11 * The above copyright notice and this permission notice shall be included in
 12 * all copies or substantial portions of the Software.
 13 *
 14 * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
 15 * IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
 16 * FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
 17 * AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
 18 * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
 19 * OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN
 20 * THE SOFTWARE.
 21 *
 22 * ************************************************************************ */
 23
 24#include "aoclsparse.h"
 25
 26#include <math.h>
 27#include <stdio.h>
 28
 29#define M 5
 30#define N 5
 31#define NNZ 8
 32
 33int main(void)
 34{
 35    aoclsparse_status     status;
 36    aoclsparse_matrix     A     = NULL;
 37    aoclsparse_mat_descr  descr = NULL;
 38    aoclsparse_operation  trans = aoclsparse_operation_none;
 39    aoclsparse_index_base base  = aoclsparse_index_base_zero;
 40
 41    double alpha = 1.0;
 42    double beta  = 0.0;
 43
 44    // Input matrix
 45    //  1  0  0  2  0
 46    //  0  3  0  0  0
 47    //  0  0  4  0  0
 48    //  0  5  0  6  7
 49    //  0  0  0  0  8
 50    aoclsparse_int csr_row_ptr[M + 1] = {0, 2, 3, 4, 7, 8};
 51    aoclsparse_int csr_col_ind[NNZ]   = {0, 3, 1, 2, 1, 3, 4, 4};
 52    double         csr_val[NNZ]       = {1, 2, 3, 4, 5, 6, 7, 8};
 53
 54    // Input vectors
 55    double x[N] = {1.0, 2.0, 3.0, 4.0, 5.0};
 56    double y[M];
 57
 58    // Expected reference result & tolerance
 59    double         y_exp[M] = {9.0, 6.0, 12.0, 69.0, 40.0};
 60    double         tol      = 1e-12;
 61    aoclsparse_int oki, ok = 1;
 62
 63    // Print aoclsparse version
 64    printf("%s\n", aoclsparse_get_version());
 65
 66    // Create matrix descriptor
 67    // aoclsparse_create_mat_descr sets aoclsparse_matrix_type to aoclsparse_matrix_type_general
 68    // and aoclsparse_index_base to aoclsparse_index_base_zero.
 69    aoclsparse_create_mat_descr(&descr);
 70
 71    // Initialise sparse matrix
 72    status = aoclsparse_create_dcsr(&A, base, M, N, NNZ, csr_row_ptr, csr_col_ind, csr_val);
 73    if(status != aoclsparse_status_success)
 74    {
 75        printf("Error while creating a sparse matrix, status = %i\n", status);
 76        return 1;
 77    }
 78
 79    // Hint the system what operation to expect // to identify hint id(which routine is to be executed, destroyed later)
 80    status = aoclsparse_set_mv_hint(A, trans, descr, 1);
 81    if(status != aoclsparse_status_success)
 82    {
 83        printf("Error while hinting operation, status = %i\n", status);
 84        return 1;
 85    }
 86
 87    // Optimize the matrix
 88    status = aoclsparse_optimize(A);
 89    if(status != aoclsparse_status_success)
 90    {
 91        printf("Error while optimizing the matrix, status = %i\n", status);
 92        return 1;
 93    }
 94
 95    // Invoke SPMV API (double precision)
 96    printf("Invoking aoclsparse_dmv...\n");
 97    status = aoclsparse_dmv(trans, &alpha, A, descr, x, &beta, y);
 98    if(status != aoclsparse_status_success)
 99    {
100        printf("Error while computing SPMV, status = %i\n", status);
101        return 1;
102    }
103
104    // Print and check the results
105    printf("Output vector:\n");
106    for(aoclsparse_int i = 0; i < M; i++)
107    {
108        oki = fabs(y[i] - y_exp[i]) <= tol;
109        printf("  %lf %c\n", y[i], oki ? ' ' : '!');
110        ok = ok && oki;
111    }
112
113    aoclsparse_destroy_mat_descr(descr);
114    aoclsparse_destroy(&A);
115    return ok ? 0 : 2;
116}
Parameters:
Return values:
  • aoclsparse_status_success – The operation completed successfully.

  • aoclsparse_status_invalid_size – The value of m, n or nnz is invalid.

  • aoclsparse_status_invalid_pointerdescr, alpha, internal structures related to the sparse matrix A, x, beta or y has an invalid pointer.

  • aoclsparse_status_not_implemented – The requested functionality is not implemented.