[petsc-users] PETSC_NULL_OBJECT gets corrupt after call to MatNestGetISs in fortran

Klaij, Christiaan C.Klaij at marin.nl
Thu Feb 12 03:02:10 CST 2015


Using petsc-3.5.3, I noticed that PETSC_NULL_OBJECT gets corrupt after calling MatNestGetISs in fortran. Here's a small example:

$ cat fieldsplittry2.F90
program fieldsplittry2

  use petscksp
  implicit none
#include <finclude/petsckspdef.h>

  PetscErrorCode :: ierr
  PetscInt       :: size,i,j,start,end,n=4,numsplit=1
  PetscScalar    :: zero=0.0,one=1.0
  Vec            :: diag3,x,b
  Mat            :: A,subA(4),myS
  PC             :: pc,subpc(2)
  KSP            :: ksp,subksp(2)
  IS             :: isg(2)

  call PetscInitialize(PETSC_NULL_CHARACTER,ierr); CHKERRQ(ierr)
  call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr); CHKERRQ(ierr);

  ! vectors
  call VecCreateMPI(MPI_COMM_WORLD,3*n,PETSC_DECIDE,diag3,ierr); CHKERRQ(ierr)
  call VecSet(diag3,one,ierr); CHKERRQ(ierr)

  call VecCreateMPI(MPI_COMM_WORLD,4*n,PETSC_DECIDE,x,ierr); CHKERRQ(ierr)
  call VecSet(x,zero,ierr); CHKERRQ(ierr)

  call VecDuplicate(x,b,ierr); CHKERRQ(ierr)
  call VecSet(b,one,ierr); CHKERRQ(ierr)

  ! matrix a00
  call MatCreateAIJ(MPI_COMM_WORLD,3*n,3*n,PETSC_DECIDE,PETSC_DECIDE,1,PETSC_NULL_INTEGER,0,PETSC_NULL_INTEGER,subA(1),ierr);CHKERRQ(ierr)
  call MatDiagonalSet(subA(1),diag3,INSERT_VALUES,ierr);CHKERRQ(ierr)
  call MatAssemblyBegin(subA(1),MAT_FINAL_ASSEMBLY,ierr);CHKERRQ(ierr)
  call MatAssemblyEnd(subA(1),MAT_FINAL_ASSEMBLY,ierr);CHKERRQ(ierr)

  ! matrix a01
  call MatCreateAIJ(MPI_COMM_WORLD,3*n,n,PETSC_DECIDE,PETSC_DECIDE,1,PETSC_NULL_INTEGER,1,PETSC_NULL_INTEGER,subA(2),ierr);CHKERRQ(ierr)
  call MatGetOwnershipRange(subA(2),start,end,ierr);CHKERRQ(ierr);
  do i=start,end-1
     j=mod(i,size*n)
     call MatSetValue(subA(2),i,j,one,INSERT_VALUES,ierr);CHKERRQ(ierr)
  end do
  call MatAssemblyBegin(subA(2),MAT_FINAL_ASSEMBLY,ierr);CHKERRQ(ierr)
  call MatAssemblyEnd(subA(2),MAT_FINAL_ASSEMBLY,ierr);CHKERRQ(ierr)

  ! matrix a10
  call  MatTranspose(subA(2),MAT_INITIAL_MATRIX,subA(3),ierr);CHKERRQ(ierr)

  ! matrix a11 (empty)
  call MatCreateAIJ(MPI_COMM_WORLD,n,n,PETSC_DECIDE,PETSC_DECIDE,0,PETSC_NULL_INTEGER,0,PETSC_NULL_INTEGER,subA(4),ierr);CHKERRQ(ierr)
  call MatAssemblyBegin(subA(4),MAT_FINAL_ASSEMBLY,ierr);CHKERRQ(ierr)
  call MatAssemblyEnd(subA(4),MAT_FINAL_ASSEMBLY,ierr);CHKERRQ(ierr)

  ! nested mat [a00,a01;a10,a11]
  call MatCreateNest(MPI_COMM_WORLD,2,PETSC_NULL_OBJECT,2,PETSC_NULL_OBJECT,subA,A,ierr);CHKERRQ(ierr)
  call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr);CHKERRQ(ierr)
  call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr);CHKERRQ(ierr)
print *, PETSC_NULL_OBJECT
  call MatNestGetISs(A,isg,PETSC_NULL_OBJECT,ierr);CHKERRQ(ierr);
print *, PETSC_NULL_OBJECT

  call PetscFinalize(ierr)

end program fieldsplittry2
$ ./fieldsplittry2
                     0
              39367824
$


dr. ir. Christiaan Klaij
CFD Researcher
Research & Development
E mailto:C.Klaij at marin.nl
T +31 317 49 33 44


MARIN
2, Haagsteeg, P.O. Box 28, 6700 AA Wageningen, The Netherlands
T +31 317 49 39 11, F +31 317 49 32 45, I www.marin.nl



More information about the petsc-users mailing list