[petsc-users] column index in MatSetValues()

Barry Smith bsmith at mcs.anl.gov
Wed Dec 1 15:48:47 CST 2010


  Humm, the problem is still very likely related to a miss-match between 4 byte and 8 byte integers.

   You should just install PETSc yourself (then you have control over it, giving control to someone else whenever doing scientific computing is always dangerous).

   Installing PETSc is usually no big deal. If you have problems send configure.log and make.log to petsc-maint at mcs.anl.gov


   Barry

On Dec 1, 2010, at 3:44 PM, Peter Wang wrote:

> Thanks,
>  
>  I changed the '1' to PetscInt ione, However, the error still comes out.
>  
>  do II=Istart,Iend-1
>         mone=II+1 !(Coef,snr,rnr are 1-based row and column numbers, shifting them to 0-based)
>             write(*,'(1a,4i7)')'II=',II,mone,snr(mone),rnr(mone)
>         call MatSetValues(A_Petsc,ione,snr(mone),ione,rnr(mone),Coef(mone),INSERT_VALUES,ierr)  ! PetscInt  ione and mone; PetscInt snr(n_nz),rnr(n_nz) PetscReal Coef(n_nz)
>                                                ^^         ^^^
>     enddo 
> 
>  
> BTW, I am running the code on the clusters of supurcomputer. Where the option  ' --with-64-bit-indices=1' shold I find and remove? 
>  
> !===The modified code is ==
> program Debug_PETSc_MatCreate_20101130
> implicit none
> !
> #include "finclude/petscsys.h"
> #include "finclude/petscvec.h"
> #include "finclude/petscmat.h"
> #include "finclude/petscpc.h"
> #include "finclude/petscksp.h"
> ! Variables
> !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
> !  PETSc Variables
> !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
>       real*8  norm
>       PetscInt  i,j,II,JJ,its   !,m,n
>       PetscInt  Istart,Iend,ione,mone
>       PetscErrorCode ierr
>       PetscMPIInt     myid,numprocs
>       PetscTruth  flg
>       PetscScalar v,one,neg_one
>       Vec         x,b,u
>       Mat         A_petsc 
>       KSP         ksp
>       PetscInt,parameter::n_nz=4
> !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
> !  Other Variables
> !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
> PetscInt::snr(n_nz),rnr(n_nz)
> !parameter::n_nz=4
> PetscReal::Coef(n_nz)
> data Coef /1.,   2.,    3. ,  4./
> data snr  /0,    1,     2,    3/
> data rnr  /0,    1 ,    2,    3/
> ! Body of Debug_PETSc_MatCreate_20101130
> ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
> !                 Beginning of program
> ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
>     call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
>     call MPI_Comm_rank(PETSC_COMM_WORLD,myid,ierr)
>     call MPI_Comm_size(PETSC_COMM_WORLD,numprocs,ierr)
>         write(*,"('snr=',4i4)")snr     
>         write(*,"('rnr=',4i4)")rnr     
>     call MatCreate(PETSC_COMM_WORLD,A_Petsc,ierr)
>     call MatSetSizes(A_Petsc,PETSC_DECIDE,PETSC_DECIDE,n_nz,n_nz,ierr)    !n_nz-1???
>     call MatSetFromOptions(A_Petsc,ierr)
>     !    write(*,*)A_petsc
>     call MatGetOwnershipRange(A_Petsc,Istart,Iend,ierr)
>     
>     write(*,'(1a,1i7,1a,1i7)')  &
>             '8.....Check after MatGetOwnershipRange() Istart=',Istart,'  Iend=',Iend
>     do II=Istart,Iend-1
>         mone=II+1 !(Coef,snr,rnr are 1-based row and column numbers, shifting them to 0-based)
>             write(*,'(1a,4i7)')'II=',II,mone,snr(mone),rnr(mone)
>         call MatSetValues(A_Petsc,ione,snr(mone),ione,rnr(mone),Coef(mone),INSERT_VALUES,ierr)
>     enddo 
>     
>         write(*,'(1a)')'9.....Check after MatSetValues()'
>     call MatAssemblyBegin(A_petsc,MAT_FINAL_ASSEMBLY,ierr)
>     call MatAssemblyEnd(A_Petsc,MAT_FINAL_ASSEMBLY,ierr)
>         write(*,'(1a)')'10.....Check after MatCreate()'
>     call MatView(A_Petsc,PETSC_VIEWER_STDOUT_WORLD,ierr)
> !      call KSPDestroy(ksp,ierr)
> !      call VecDestroy(u,ierr)
> !      call VecDestroy(x,ierr)
> !      call VecDestroy(b,ierr)
>     call MatDestroy(A_petsc,ierr)
>     call PetscFinalize(ierr)
> end program Debug_PETSc_MatCreate_20101130
> 
>  
> > From: bsmith at mcs.anl.gov
> > Date: Wed, 1 Dec 2010 08:06:19 -0600
> > To: petsc-users at mcs.anl.gov
> > Subject: Re: [petsc-users] column index in MatSetValues()
> > 
> > 
> > 
> > --with-64-bit-indices=1
> > 
> > You only need this option if you are solving problems with over 2 billion unknowns! I recommend removing it otherwise, it wastes memory and slows performance slightly.
> > 
> > > MatSetValues(A_Petsc,1,snr(Ione),1,rnr(Ione),Coef(Ione)
> > ^^^^ ^^^^^
> > 
> > --with-64-bit-indices means ALL integers passed to PETSc MUST be 64 bit, but here you are passing the integer 1 as a "regular" 32 bit integer. You need to declare it as a PetscInt, for example
> > 
> > PetscInt mone
> > mone = 1
> > > MatSetValues(A_Petsc,mone,snr(Ione),mone,rnr(Ione),Coef(Ione)
> > 
> > but better just build PETSc without the --with-64-bit-indices
> > 
> > Barry
> > 
> > 



More information about the petsc-users mailing list