[petsc-users] column index in MatSetValues()

Peter Wang pengxwang at hotmail.com
Thu Dec 2 12:17:37 CST 2010


  I sent two emails for replying thie topic. However, I didn't get the email of myself from petsc-users-bounces at mcs.anl.gov .  I am wondering if the email system has something wrong? Sorry if the resent email bohters anyone.
 
  In the new version of the code I defined PetscInt II,JJ,ione, mone,snr[] and rnr[], PetscReal coef[], and modified the following portion. However, the error is still there. Is there any ohter reason I didn't figure out?
 
 
 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 loop=======
    
  do I=Istart,Iend-1
        mone=I+1 !(Coef,snr,rnr are 1-based row and column numbers, shifting them to 0-based)
        II=snr(mone)
        JJ=rnr(mone)
        v=coef(mone)
            write(*,'(1a,4i7)')'II=',II,mone,snr(mone),rnr(mone)
        call MatSetValues(A_Petsc,ione,II,ione,JJ,v,INSERT_VALUES,ierr)
    enddo 

! ============the whole program code modified====================
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 I=Istart,Iend-1
        mone=I+1 !(Coef,snr,rnr are 1-based row and column numbers, shifting them to 0-based)
        II=snr(mone)
        JJ=rnr(mone)
        v=coef(mone)
            write(*,'(1a,4i7)')'II=',II,mone,snr(mone),rnr(mone)
        call MatSetValues(A_Petsc,ione,II,ione,JJ,v,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
 
!===================End of the code============================

 

 

> 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
> 
> 
> On Nov 30, 2010, at 10:33 PM, Peter Wang wrote:
> 
> > I am trying to create a matrix and insert values to it. The martix is supposed to be as following:
> > 
> > 1 0 0 0
> > 0 2 0 0
> > 0 0 3 0
> > 0 0 0 4
> > 
> > array coef[] is the diagonal value of the matrix, 
> > snr[] is the index of the row, rnr[] is the index of column.
> > 
> > However, I always get the wrong results. It shows the Column too large: col 4607182418800017408 max 3! I cheked the value of rnr[]. The output snr and rnr is correct:
> > snr= 0 1 2 3
> > rnr= 0 1 2 3
> > 
> > It seems there is something wrong when MatSetValues() is called. Following is a part of the error information. The information is shown at each loop of do II=Istart,Iend-1
> > 
> > The output (if any) follows:
> > snr= 0 1 2 3
> > rnr= 0 1 2 3
> > 8.....Check after MatGetOwnershipRange() Istart= 0 Iend= 4
> > II= 0 1 0 0
> > [0]PETSC ERROR: --------------------- Error Message ------------------------------------
> > [0]PETSC ERROR: Argument out of range!
> > [0]PETSC ERROR: Column too large: col 4607182418800017408 max 3!
> > [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: Debug_PETSc_MatCreate_20101130 on a linux-gnu named compute-1-35.hpc.local.uwm by pwang_a Tue Nov 30 22:27:03 2010
> > [0]PETSC ERROR: Libraries linked from /sharedapps/uwm/ceas/gcc-4.4.3/petsc/3.1-p5-v1/lib
> > [0]PETSC ERROR: Configure run at Fri Oct 8 12:59:16 2010
> > [0]PETSC ERROR: Configure options --prefix=/sharedapps/uwm/ceas/gcc-4.4.3/petsc/3.1-p5-v1 --with-mpi-dir=/sharedapps/uwm/common/gcc-4.4.3/openmpi/1.3.2-v1 --with-blas-lapack-dir=/sharedapps/uwm/ceas/gcc-4.4.3/lapack/3.2.2-v1/lib --with-64-bit-indices=1 --with-64-bit-pointers=1 --with-large-file-io=1 --with-x=0
> > [0]PETSC ERROR: ------------------------------------------------------------------------
> > [0]PETSC ERROR: MatSetValues_SeqAIJ() line 193 in src/mat/impls/aij/seq/aij.c
> > [0]PETSC ERROR: MatSetValues() line 992 in src/mat/interface/matrix.c
> > 
> > 
> > 
> > !The code is as following:
> > !=============================
> > 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
> > 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
> > !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
> > !parameter::n_nz=4
> > Real*8::Coef(n_nz)
> > PetscInt::snr(n_nz),rnr(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
> > ione=II+1 !(Coef,snr,rnr are 1-based row and column numbers, shifting them to 0-based)
> > write(*,'(1a,4i7)')'II=',II,ione,snr(ione),rnr(ione) !output snr and rnr for error check
> > call MatSetValues(A_Petsc,1,snr(Ione),1,rnr(Ione),Coef(Ione),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
> > !=====================================
> 
 		 	   		  
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20101202/581ad5af/attachment-0001.htm>


More information about the petsc-users mailing list