[petsc-users] From CSR format to MatSetValues

Matthew Knepley knepley at gmail.com
Sun Jun 10 07:49:13 CDT 2012


On Sun, Jun 10, 2012 at 4:48 PM, w_ang_temp <w_ang_temp at 163.com> wrote:

> Hello, Barry
>
>     I try to use MatCreateSeqAIJWithArrays to deal with my problem. Since
> there is no example in the doc, I find some examples in the web and use it.
> And I get some errors.
>
>     When the processor is 1(mpiexec -n 1 ./ex4f), the results are ok. While
> it is more than one, I get the errors. So I think I miss some piece of
> information about MatCreateSeqAIJWithArrays.
>

1) Check the error codes. Always!


>     Thanks.
>                                           Jim
>
> ----------------------------code-------------------------------------------
>       call MatCreate(PETSC_COMM_WORLD,A,ierr)
>       call MatSetType(A,MATAIJ,ierr)
>       call MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n,ierr)
>       call MatSetFromOptions(A,ierr)
>       call MatGetOwnershipRange(A,Istart,Iend,ierr)
>


>       !NROWIN,NNZIJ,values:the CSR variables of the main function
>       !NROWIN,NNZIJ:one-based,so get zero-based variables rowIndex,columns
>       allocate(rowIndex(JFLNG+1))
>       allocate(columns(MAXNNZ))
>       rowIndex=NROWIN-1
>       columns=NNZIJ-1
>
>       call MatCreateSeqAIJWithArrays(MPI_COMM_WORLD,m,n,rowIndex,columns,
>      &       values,A, ierr)
>

2) This can only be called in serial. In parallel, you need
MatCreateMPIAIJWithArrays()

3) I still believe you are really doing the wrong thing. The right solution
is to set one row at
    a time. I will bet dollars to doughnuts that the assembly time is not
noticeable in your program.

     Matt


>
>       call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
>       call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)
> ----------------------------code-------------------------------------------
>
> ----------------------------Error Message----------------------------------
> .........
> [1]PETSC ERROR: --------------------- Error Message
> ------------------------------------
> [1]PETSC ERROR: Object is in wrong state!
> [1]PETSC ERROR: Mat object's type is not set: Argument # 1!
> [1]PETSC ERROR:
> ------------------------------------------------------------------------
> [1]PETSC ERROR: Petsc Release Version 3.2.0, Patch 7, Thu Mar 15 09:30:51
> CDT 2012
> [1]PETSC ERROR: See docs/changes/index.html for recent updates.
> [1]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
> [1]PETSC ERROR: See docs/index.html for manual pages.
> [1]PETSC ERROR:
> ------------------------------------------------------------------------
> [1]PETSC ERROR: ./ex4f on a arch-linu named ubuntu by geo Sun Jun 10
> 01:28:48 2012
> [1]PETSC ERROR: Libraries linked from
> /home/geo/soft/petsc/petsc-3.2-p7/arch-linux2-c-debug/lib
> [1]PETSC ERROR: Configure run at Sat Jun  2 00:53:03 2012
> [1]PETSC ERROR: Configure options --with-mpi-dir=/home/geo/soft/mpich2
> --download-f-blas-lapack=1 --with-x=1
> [1]PETSC ERROR:
> ------------------------------------------------------------------------
> [1]PETSC ERROR: MatScaleSystem() line 920 in src/mat/interface/matrix.c
> [1]PETSC ERROR: PCPreSolve() line 1376 in src/ksp/pc/interface/precon.c
> [1]PETSC ERROR: KSPSolve() line 404 in src/ksp/ksp/interface/itfunc.c
> [1]PETSC ERROR: --------------------- Error Message
> ------------------------------------
> [1]PETSC ERROR: Argument out of range!
> [1]PETSC ERROR: Comm must be of size 1!
> [1]PETSC ERROR:
> ------------------------------------------------------------------------
> [1]PETSC ERROR: Petsc Release Version 3.2.0, Patch 7, Thu Mar 15 09:30:51
> CDT 2012
> [1]PETSC ERROR: See docs/changes/index.html for recent updates.
> [1]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
> [1]PETSC ERROR: See docs/index.html for manual pages.
> [1]PETSC ERROR:
> ------------------------------------------------------------------------
> [1]PETSC ERROR: ./ex4f on a arch-linu named ubuntu by geo Sun Jun 10
> 01:28:48 2012
> [1]PETSC ERROR: Libraries linked from
> /home/geo/soft/petsc/petsc-3.2-p7/arch-linux2-c-debug/lib
> [1]PETSC ERROR: Configure run at Sat Jun  2 00:53:03 2012
> [1]PETSC ERROR: Configure options --with-mpi-dir=/home/geo/soft/mpich2
> --download-f-blas-lapack=1 --with-x=1
> [1]PETSC ERROR:
> ------------------------------------------------------------------------
> [1]PETSC ERROR: MatCreate_SeqAIJ() line 3794 in src/mat/impls/aij/seq/aij.c
> [1]PETSC ERROR: MatSetType() line 73 in src/mat/interface/matreg.c
> [1]PETSC ERROR: MatCreateSeqAIJWithArrays() line 4171 in
> src/mat/impls/aij/seq/aij.c
> [1]PETSC ERROR: --------------------- Error Message
> ------------------------------------
> [1]PETSC ERROR: Object is in wrong state!
> [1]PETSC ERROR: Mat object's type is not set: Argument # 1!
> [1]PETSC ERROR:
> ------------------------------------------------------------------------
> [1]PETSC ERROR: Petsc Release Version 3.2.0, Patch 7, Thu Mar 15 09:30:51
> CDT 2012
> [1]PETSC ERROR: See docs/changes/index.html for recent updates.
> [1]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
> [1]PETSC ERROR: See docs/index.html for manual pages.
> [1]PETSC ERROR:
> ------------------------------------------------------------------------
> [1]PETSC ERROR: ./ex4f on a arch-linu named ubuntu by geo Sun Jun 10
> 01:28:48 2012
> [1]PETSC ERROR: Libraries linked from
> /home/geo/soft/petsc/petsc-3.2-p7/arch-linux2-c-debug/lib
> [1]PETSC ERROR: Configure run at Sat Jun  2 00:53:03 2012
> [1]PETSC ERROR: Configure options --with-mpi-dir=/home/geo/soft/mpich2
> --download-f-blas-lapack=1 --with-x=1
> [1]PETSC ERROR:
> ------------------------------------------------------------------------
> [1]PETSC ERROR: MatAssemblyBegin() line 4760 in src/mat/interface/matrix.c
> [1]PETSC ERROR: --------------------- Error Message
> ------------------------------------
> [1]PETSC ERROR: Object is in wrong state!
> [1]PETSC ERROR: Mat object's type is not set: Argument # 1!
> [1]PETSC ERROR:
> ------------------------------------------------------------------------
> [1]PETSC ERROR: Petsc Release Version 3.2.0, Patch 7, Thu Mar 15 09:30:51
> CDT 2012
> [1]PETSC ERROR: See docs/changes/index.html for recent updates.
> [1]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
> [1]PETSC ERROR: See docs/index.html for manual pages.
> [1]PETSC ERROR:
> ------------------------------------------------------------------------
> [1]PETSC ERROR: ./ex4f on a arch-linu named ubuntu by geo Sun Jun 10
> 01:28:48 2012
> [1]PETSC ERROR: Libraries linked from
> /home/geo/soft/petsc/petsc-3.2-p7/arch-linux2-c-debug/lib
> [1]PETSC ERROR: Configure run at Sat Jun  2 00:53:03 2012
> [1]PETSC ERROR: Configure options --with-mpi-dir=/home/geo/soft/mpich2
> --download-f-blas-lapack=1 --with-x=1
> [1]PETSC ERROR:
> ------------------------------------------------------------------------
> [1]PETSC ERROR: MatAssemblyEnd() line 4937 in src/mat/interface/matrix.c
> [1]PETSC ERROR: --------------------- Error Message
> ------------------------------------
> [1]PETSC ERROR: Object is in wrong state!
> [1]PETSC ERROR: Mat object's type is not set: Argument # 1!
> [1]PETSC ERROR:
> ------------------------------------------------------------------------
> [1]PETSC ERROR: Petsc Release Version 3.2.0, Patch 7, Thu Mar 15 09:30:51
> CDT 2012
> [1]PETSC ERROR: See docs/changes/index.html for recent updates.
> [1]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
> [1]PETSC ERROR: See docs/index.html for manual pages.
> [1]PETSC ERROR:
> ------------------------------------------------------------------------
> [1]PETSC ERROR: ./ex4f on a arch-linu named ubuntu by geo Sun Jun 10
> 01:28:48 2012
> [1]PETSC ERROR: Libraries linked from
> /home/geo/soft/petsc/petsc-3.2-p7/arch-linux2-c-debug/lib
> [1]PETSC ERROR: Configure run at Sat Jun  2 00:53:03 2012
> [1]PETSC ERROR: Configure options --with-mpi-dir=/home/geo/soft/mpich2
> --download-f-blas-lapack=1 --with-x=1
> [1]PETSC ERROR:
> ------------------------------------------------------------------------
> [1]PETSC ERROR: MatScaleSystem() line 920 in src/mat/interface/matrix.c
> [1]PETSC ERROR: PCPreSolve() line 1376 in src/ksp/pc/interface/precon.c
> [1]PETSC ERROR: KSPSolve() line 404 in src/ksp/ksp/interface/itfunc.c
> ----------------------------Error Message----------------------------------
>
>
>
> At 2012-06-09 06:10:50,"Barry Smith" <bsmith at mcs.anl.gov> wrote:
> >
> >http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatCreateSeqAIJWithArrays.html
> >http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatCreateMPIAIJWithArrays.html
> >
> >
> >
> >On Jun 8, 2012, at 1:31 AM, w_ang_temp wrote:
> >
> >> Hello,
> >>     When solving a linear system Ax=b, the first step is assigning values to
> >> matrix A. In my program, the subroutine PETSCSOLVE, which is used to slove Ax=b
> >> with PETSc, gets the CSR format matrix(values, columns, rowIndex) to set values
> >> to PETSc Mat A.
> >>     The variables 'values'、'columns'、'rowIndex' belong to the main function.
> >> They are used to represent a matrix in CSR format. The following table describes
> >> the arrays in terms of the values, row, and column positions of the non-zero
> >> elements in a sparse matrix.
> >>     values: A real or complex array that contains the non-zero elements of a
> >>             sparse matrix. The non-zero elements are mapped into the values
> >>             array using the row-major upper triangular storage mapping described
> >>             above.
> >>     columns: Element i of the integer array columns is the number of the column
> >>              that contains the i-th element in the values array.
> >>     rowIndex: Element j of the integer array rowIndex gives the index of the
> >>               element in the values array that is first non-zero element in a row j.
> >>
> >> codes:
> >> ---------Set Values to A From CSR Format(values,rowIndex,columns)-------------------
> >>       !values:Non-Symmetric Matrix
> >>       ione = 1
> >>       do 10,II=Istart+1,Iend
> >>         !non-zero numbers of this row
> >>         rowNum=rowIndex(II+1)-rowIndex(II)
> >>         do 11,JJ=1,rowNum
> >>
> >>             !elemnt index of the values/columns
> >>             kValNn=rowIndex(II)+JJ-1
> >>             !column index of this element
> >>             nCol=columns(kValNn)-1
> >>
> >>             !Setdata:II-RowIndex,nCol-ColIndex
> >>             nRow=II-1
> >>             call MatSetValues(A,ione,nRow,ione,nCol,values(kValNn),INSERT_VALUES,ierr)
> >>   11    continue
> >>   10  continue
> >> --------------------------------------------------------------------------------------
> >>
> >>     As we can see, it has to loop a number of times because it can only set one data per
> >> time. And this leads to low efficiency.
> >>     So is there a better way to do this?
> >>     Thank you.
> >>                                       Jim
> >>
> >>
> >>
> >
>
>
>
>


-- 
What most experimenters take for granted before they begin their
experiments is infinitely more interesting than any results to which their
experiments lead.
-- Norbert Wiener
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20120610/1b91a526/attachment-0001.html>


More information about the petsc-users mailing list