[petsc-users] problems with MatCreateSeqAIJWithArrays

Nisoli Isaia isaia.nisoli at gmail.com
Mon Nov 28 19:03:26 CST 2011


I have a small program that takes as input a file in UBLAS serialization
format of a UBLAS compressed matrix and outputs a Petsc binary file, using
MatCreateSeqAIJWithArrays.

int convertPETSc(srmatrix &sA, Mat * const Aaddr)
{
    PetscErrorCode ierr;
    PetscInt size=(sA.A).size1();
    unsigned int totnnz=(sA.A).nnz();

    (sA.A).complete_index1_data();
    long unsigned int *row_ptr =(sA.A).index1_data().begin();
    long unsigned int *col_ptr =(sA.A).index2_data().begin();
    double * value_ptr = (sA.A).value_data().begin();
    unsigned int sizerow_ptr=(sA.A).index1_data().size();
    unsigned int sizecol_ptr=(sA.A).index2_data().size();
    PetscInt* rowindices;
    PetscInt* colindices;
    PetscScalar* values;
    rowindices=new PetscInt[sizerow_ptr];
    colindices=new PetscInt[sizecol_ptr];
    values=new PetscScalar[totnnz];
    for (unsigned int i=0;i<sizerow_ptr;i++)
rowindices[i]=PetscInt(row_ptr[i]);
    for (unsigned int i=0;i<sizecol_ptr;i++)
colindices[i]=PetscInt(col_ptr[i]);
    for (unsigned int i=0;i<totnnz;i++) values[i]=PetscScalar(value_ptr[i]);

ierr=MatCreateSeqAIJWithArrays(PETSC_COMM_WORLD,size+1,size+1,rowindices,colindices,values,Aaddr);CHKERRQ(ierr);
    ierr=MatAssemblyBegin(*Aaddr,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
    ierr=MatAssemblyEnd(*Aaddr,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
    return 0;
}

It writes the output without any problem; the problem is  that when I load
the file I get the following error message.

A strange issue is that, if I transpose the matrix before outputting it, I
have no problem.

Please note that nnz for the example matrix is 2227712; I don't understand
why it thinks it is -2227712.

Thanks in advance
Isaia

[0]PETSC ERROR: --------------------- Error Message
------------------------------------
[0]PETSC ERROR: Argument out of range!
[0]PETSC ERROR: nnz cannot be less than 0: local row 524288 value -2227712!
[0]PETSC ERROR:
------------------------------------------------------------------------
[0]PETSC ERROR: Petsc Release Version 3.1.0, Patch 5, Mon Sep 27 11:51:54
CDT 2010
[0]PETSC ERROR: See docs/changes/index.html for recent updates.
[0]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
[0]PETSC ERROR: See docs/index.html for manual pages.
[0]PETSC ERROR:
------------------------------------------------------------------------
[0]PETSC ERROR: ./convertPetsc on a linux-gnu named aiace2 by isaia Tue Nov
29 01:49:57 2011
[0]PETSC ERROR: Libraries linked from
/build/buildd/petsc-3.1.dfsg/linux-gnu-c-opt/lib
[0]PETSC ERROR: Configure run at Mon Jan  3 21:43:51 2011
[0]PETSC ERROR: Configure options --with-shared --with-debugging=0
--useThreads 0 --with-clanguage=C++ --with-c-support
--with-fortran-interfaces=1 --with-mpi-dir=/usr/lib/openmpi
--with-mpi-shared=1 --with-blas-lib=-lblas --with-lapack-lib=-llapack
--with-blacs=1 --with-blacs-include=/usr/include
--with-blacs-lib="[/usr/lib/libblacsCinit-openmpi.so,/usr/lib/libblacs-openmpi.so]"
--with-scalapack=1 --with-scalapack-include=/usr/include
--with-scalapack-lib=/usr/lib/libscalapack-openmpi.so --with-mumps=1
--with-mumps-include=/usr/include
--with-mumps-lib="[/usr/lib/libdmumps.so,/usr/lib/libzmumps.so,/usr/lib/libsmumps.so,/usr/lib/libcmumps.so,/usr/lib/libmumps_common.so,/usr/lib/libpord.so]"
--with-umfpack=1 --with-umfpack-include=/usr/include/suitesparse
--with-umfpack-lib="[/usr/lib/libumfpack.so,/usr/lib/libamd.so]"
--with-spooles=1 --with-spooles-include=/usr/include/spooles
--with-spooles-lib=/usr/lib/libspooles.so --with-hypre=1
--with-hypre-dir=/usr --with-scotch=1
--with-scotch-include=/usr/include/scotch
--with-scotch-lib=/usr/lib/libscotch.so --with-hdf5=1 --with-hdf5-dir=/usr
[0]PETSC ERROR:
------------------------------------------------------------------------
[0]PETSC ERROR: MatSeqAIJSetPreallocation_SeqAIJ() line 3002 in
src/mat/impls/aij/seq/aij.c
[0]PETSC ERROR: MatLoad_SeqAIJ() line 3566 in src/mat/impls/aij/seq/aij.c
[0]PETSC ERROR: MatLoad() line 124 in src/mat/utils/matio.c
[0]PETSC ERROR: User provided function() line 119 in convertPetsc.cxx
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20111129/5e1bc158/attachment.htm>


More information about the petsc-users mailing list