Intel® oneAPI Math Kernel Library
Ask questions and share information with other developers who use Intel® Math Kernel Library.

mkl_sparse_z_export_csr of Inspector-executor Sparse BLAS return bizarre results.

Dan_Ghiocel
Beginner
367 Views

Dear MKL experts,

       I'm testing the function,mkl_sparse_z_export_csr, in Inspector-executor Sparse BLAS, to convert COO format into CSR format and output the CSR format. The function returns BIZARRE results as shown in the below. I guess I did something wrong. I attach the simple code at the end. Please tell what the correct way to do it.

Dan

!   Sparse representation of the matrix A
!
!                 |   1       -1      0   -3     0   |
!                 |  -2        5      0    0     0   |
!   A    =     |   0        0      4    6     4    |,
!                 |  -4        0      2    7     0   |
!                 |   0        8      0    0    -5   |
!

   info = mkl_sparse_z_export_csr (csrA, i, nrow, ncol, rows_start,rows_end,col_indx,csrA_value)

 OUTPUT DATA FOR in CSR Format =
 rows_start(1:M) =      8546432  1919251285  1631870067  1883331694  1952531568
 rows_end(1:M)  =      8546436  1952531568  1867275361  1550606691  1886217556
 col_indx(1:NNZ) =     8546688  1279350084  1313817944  1398079538  1095651909  1144866125
                                 1426091617  1347568979  1229344594  1128088908...














 

 

0 Kudos
1 Solution
MariaZh
Employee
367 Views

Hello Dan,

Here is the example which should help to resolve incorrect behavior of mkl_sparse_z_export_csr by adding mkl_sparse_z_export_csr_cf interface which is using C_PTR.
Please let me know what you find out!

Best regards,
Maria

 

PROGRAM EXPORT_CSR

    USE MKL_SPBLAS
    USE ISO_C_BINDING
    IMPLICIT NONE
!   *****************************************************************************
!   Declare new interface using C_PTR:
!   *****************************************************************************
    INTERFACE
        FUNCTION MKL_SPARSE_Z_EXPORT_CSR_CF(source,indexing,rows,cols,rows_start,rows_end,col_indx,values) &
                 BIND(C, name='MKL_SPARSE_Z_EXPORT_CSR')
            USE, INTRINSIC :: ISO_C_BINDING , ONLY : C_INT, C_DOUBLE_COMPLEX, C_PTR
            IMPORT SPARSE_MATRIX_T
            TYPE(SPARSE_MATRIX_T) , INTENT(IN) :: source
            INTEGER(C_INT), INTENT(INOUT) :: indexing
            INTEGER       , INTENT(INOUT) :: rows
            INTEGER       , INTENT(INOUT) :: cols
            TYPE(C_PTR)   , INTENT(INOUT) :: rows_start
            TYPE(C_PTR)   , INTENT(INOUT) :: rows_end
            TYPE(C_PTR)   , INTENT(INOUT) :: col_indx
            TYPE(C_PTR)   , INTENT(INOUT) :: values
            INTEGER(C_INT) MKL_SPARSE_Z_EXPORT_CSR_CF
        END FUNCTION
    END INTERFACE

!   *****************************************************************************
!   Sparse representation of the matrices A:
!   *****************************************************************************
    INTEGER   , ALLOCATABLE :: csrColInd(:), csrRowPtr(:)
    COMPLEX*16, ALLOCATABLE :: csrVal(:)
!   CSR matrix structure
    TYPE(SPARSE_MATRIX_T) csrA
!   Variables used for exporting sparse matrix
    INTEGER        :: nrows, ncols
    INTEGER(C_INT) :: indexing
    TYPE(C_PTR)    :: rows_start_c, rows_end_c, col_indx_c, values_c
    INTEGER   , POINTER :: rows_start_f(:), rows_end_f(:), col_indx_f(:)
    COMPLEX*16, POINTER :: values_f(:)
!   *****************************************************************************
!   Declaration of local variables:
!   *****************************************************************************
    INTEGER M, N, NNZ, i, j, info
    M = 5
    NNZ = 13
    ALLOCATE(csrColInd(NNZ))
    ALLOCATE(csrRowPtr(M+1))
    ALLOCATE(csrVal(NNZ))
    csrVal = (/ (1.0, 0.0), (-1.0, 0.0), (-3.0, 0.0), (-2.0, 0.0), (5.0, 0.0), (4.0, 0.0), &
                (6.0, 0.0), (4.0, 0.0), (-4.0, 0.0), (2.0, 0.0), (7.0, 0.0), (8.0, 0.0), (-5.0, 0.0) /)
    csrColInd = (/ 0,1,3,0,1,2,3,4,0,2,3,1,4 /)
    csrRowPtr = (/ 0,    3,  5,    8,   11,  13 /)

    print*,'---------------------------------------------------'
    print*,'Input matrix A:'
    do i = 1, M
        print 100,'row #',i
        do j = csrRowPtr(i)+1, csrRowPtr(i+1)
            print*,csrColInd(j),csrVal(j)
        enddo
    enddo

!   Create CSR matrix
    info = MKL_SPARSE_Z_CREATE_CSR(csrA,SPARSE_INDEX_BASE_ZERO,M,M,csrRowPtr(1),csrRowPtr(2),csrColInd,csrVal)

!   Export CSR matrix
    info = MKL_SPARSE_Z_EXPORT_CSR_CF(csrA, indexing, nrows, ncols, rows_start_c, rows_end_c, col_indx_c, values_c)

!   Converting C into Fortran pointers
    call C_F_POINTER(rows_start_c, rows_start_f, [nrows])
    call C_F_POINTER(rows_end_c  , rows_end_f  , [nrows])
    call C_F_POINTER(col_indx_c  , col_indx_f  , [rows_end_f(nrows)])
    call C_F_POINTER(values_c    , values_f    , [rows_end_f(nrows)])

    print*,'---------------------------------------------------'
    print *,'Output matrix A:'
    do i = 1, nrows
        print 100,'row #',i
        do j = rows_start_f(i)+1, rows_end_f(i)
            print*,col_indx_f(j),values_f(j)
        enddo
    enddo

!   Release internal representation of CSR matrix
    info = MKL_SPARSE_DESTROY(csrA)

    print*,'---------------------------------------------------'

100 format(A,I2)

END PROGRAM EXPORT_CSR

 

 

View solution in original post

0 Kudos
2 Replies
MariaZh
Employee
368 Views

Hello Dan,

Here is the example which should help to resolve incorrect behavior of mkl_sparse_z_export_csr by adding mkl_sparse_z_export_csr_cf interface which is using C_PTR.
Please let me know what you find out!

Best regards,
Maria

 

PROGRAM EXPORT_CSR

    USE MKL_SPBLAS
    USE ISO_C_BINDING
    IMPLICIT NONE
!   *****************************************************************************
!   Declare new interface using C_PTR:
!   *****************************************************************************
    INTERFACE
        FUNCTION MKL_SPARSE_Z_EXPORT_CSR_CF(source,indexing,rows,cols,rows_start,rows_end,col_indx,values) &
                 BIND(C, name='MKL_SPARSE_Z_EXPORT_CSR')
            USE, INTRINSIC :: ISO_C_BINDING , ONLY : C_INT, C_DOUBLE_COMPLEX, C_PTR
            IMPORT SPARSE_MATRIX_T
            TYPE(SPARSE_MATRIX_T) , INTENT(IN) :: source
            INTEGER(C_INT), INTENT(INOUT) :: indexing
            INTEGER       , INTENT(INOUT) :: rows
            INTEGER       , INTENT(INOUT) :: cols
            TYPE(C_PTR)   , INTENT(INOUT) :: rows_start
            TYPE(C_PTR)   , INTENT(INOUT) :: rows_end
            TYPE(C_PTR)   , INTENT(INOUT) :: col_indx
            TYPE(C_PTR)   , INTENT(INOUT) :: values
            INTEGER(C_INT) MKL_SPARSE_Z_EXPORT_CSR_CF
        END FUNCTION
    END INTERFACE

!   *****************************************************************************
!   Sparse representation of the matrices A:
!   *****************************************************************************
    INTEGER   , ALLOCATABLE :: csrColInd(:), csrRowPtr(:)
    COMPLEX*16, ALLOCATABLE :: csrVal(:)
!   CSR matrix structure
    TYPE(SPARSE_MATRIX_T) csrA
!   Variables used for exporting sparse matrix
    INTEGER        :: nrows, ncols
    INTEGER(C_INT) :: indexing
    TYPE(C_PTR)    :: rows_start_c, rows_end_c, col_indx_c, values_c
    INTEGER   , POINTER :: rows_start_f(:), rows_end_f(:), col_indx_f(:)
    COMPLEX*16, POINTER :: values_f(:)
!   *****************************************************************************
!   Declaration of local variables:
!   *****************************************************************************
    INTEGER M, N, NNZ, i, j, info
    M = 5
    NNZ = 13
    ALLOCATE(csrColInd(NNZ))
    ALLOCATE(csrRowPtr(M+1))
    ALLOCATE(csrVal(NNZ))
    csrVal = (/ (1.0, 0.0), (-1.0, 0.0), (-3.0, 0.0), (-2.0, 0.0), (5.0, 0.0), (4.0, 0.0), &
                (6.0, 0.0), (4.0, 0.0), (-4.0, 0.0), (2.0, 0.0), (7.0, 0.0), (8.0, 0.0), (-5.0, 0.0) /)
    csrColInd = (/ 0,1,3,0,1,2,3,4,0,2,3,1,4 /)
    csrRowPtr = (/ 0,    3,  5,    8,   11,  13 /)

    print*,'---------------------------------------------------'
    print*,'Input matrix A:'
    do i = 1, M
        print 100,'row #',i
        do j = csrRowPtr(i)+1, csrRowPtr(i+1)
            print*,csrColInd(j),csrVal(j)
        enddo
    enddo

!   Create CSR matrix
    info = MKL_SPARSE_Z_CREATE_CSR(csrA,SPARSE_INDEX_BASE_ZERO,M,M,csrRowPtr(1),csrRowPtr(2),csrColInd,csrVal)

!   Export CSR matrix
    info = MKL_SPARSE_Z_EXPORT_CSR_CF(csrA, indexing, nrows, ncols, rows_start_c, rows_end_c, col_indx_c, values_c)

!   Converting C into Fortran pointers
    call C_F_POINTER(rows_start_c, rows_start_f, [nrows])
    call C_F_POINTER(rows_end_c  , rows_end_f  , [nrows])
    call C_F_POINTER(col_indx_c  , col_indx_f  , [rows_end_f(nrows)])
    call C_F_POINTER(values_c    , values_f    , [rows_end_f(nrows)])

    print*,'---------------------------------------------------'
    print *,'Output matrix A:'
    do i = 1, nrows
        print 100,'row #',i
        do j = rows_start_f(i)+1, rows_end_f(i)
            print*,col_indx_f(j),values_f(j)
        enddo
    enddo

!   Release internal representation of CSR matrix
    info = MKL_SPARSE_DESTROY(csrA)

    print*,'---------------------------------------------------'

100 format(A,I2)

END PROGRAM EXPORT_CSR

 

 

0 Kudos
Dan_Ghiocel
Beginner
367 Views

Hi, Maria,

      Thanks a lot. It works.

  Dan

0 Kudos
Reply