Developer Reference

Contents

mkl_sparse_?_spmmd

Computes the product of two sparse matrices and stores the result as a dense matrix.

Syntax

sparse_status_t mkl_sparse_s_spmmd
(
const sparse_operation_t
operation
,
const
sparse_matrix_t
A
,
const
sparse_matrix_t
B
,
const sparse_layout_t
layout
,
float
*C
,
const MKL_INT
ldc
);
sparse_status_t mkl_sparse_d_spmmd
(
const sparse_operation_t
operation
,
const
sparse_matrix_t
A
,
const
sparse_matrix_t
B
,
const sparse_layout_t
layout
,
double
*C
,
const MKL_INT
ldc
);
sparse_status_t mkl_sparse_c_spmmd
(
const sparse_operation_t
operation
,
const
sparse_matrix_t
A
,
const
sparse_matrix_t
B
,
const sparse_layout_t
layout
,
MKL_Complex8
*C
,
const MKL_INT
ldc
);
sparse_status_t mkl_sparse_z_spmmd
(
const sparse_operation_t
operation
,
const
sparse_matrix_t
A
,
const
sparse_matrix_t
B
,
const sparse_layout_t
layout
,
MKL_Complex16
*C
,
const MKL_INT
ldc
);
Include Files
  • mkl_spblas.h
Description
The
mkl_sparse_?_spmmd
routine performs a matrix-matrix operation:
C := op(A)*B
where
A
and
B
are sparse matrices,
op
is a matrix modifier for matrix
A
, and
C
is a dense matrix.
This routine is not supported for sparse matrices in the COO format. For sparse matrices in BSR format, these combinations of (
indexing
,
block_layout
) are supported:
  • (
    SPARSE_INDEX_BASE_ZERO
    ,
    SPARSE_LAYOUT_ROW_MAJOR
    )
  • (
    SPARSE_INDEX_BASE_ONE
    ,
    SPARSE_LAYOUT_COLUMN_MAJOR
    )
Input Parameters
operation
Specifies operation
op()
on input matrix.
SPARSE_OPERATION_NON_TRANSPOSE
Non-transpose,
op(
A
) =
A
.
SPARSE_OPERATION_TRANSPOSE
Transpose,
op(
A
) =
A
T
.
SPARSE_OPERATION_CONJUGATE_TRANSPOSE
Conjugate transpose,
op(
A
) =
A
H
.
A
Handle which contains the sparse matrix
A
.
B
Handle which contains the sparse matrix
B
.
layout
Describes the storage scheme for the dense matrix:
SPARSE_LAYOUT_COLUMN_MAJOR
Storage of elements uses column major layout.
SPARSE_LAYOUT_ROW_MAJOR
Storage of elements uses row major layout.
ldC
Leading dimension of matrix
C
.
Output Parameters
C
Resulting dense matrix.
Return Values
The function returns a value indicating whether the operation was successful or not, and why.
SPARSE_STATUS_SUCCESS
The operation was successful.
SPARSE_STATUS_NOT_INITIALIZED
The routine encountered an empty handle or matrix array.
SPARSE_STATUS_ALLOC_FAILED
Internal memory allocation failed.
SPARSE_STATUS_INVALID_VALUE
The input parameters contain an invalid value.
SPARSE_STATUS_EXECUTION_FAILED
Execution failed.
SPARSE_STATUS_INTERNAL_ERROR
An error in algorithm implementation occurred.
SPARSE_STATUS_NOT_SUPPORTED
The requested operation is not supported.

Product and Performance Information

1

Performance varies by use, configuration and other factors. Learn more at www.Intel.com/PerformanceIndex.