Actual source code: ex1f90.F90

  1: !
  2: !  Description: Creates an index set based on a set of integers. Views that index set
  3: !  and then destroys it.
  4: !
  5: !

  7:       program main
  8: #include <petsc/finclude/petscis.h>
  9:       use petscis
 10:       implicit none

 12:       PetscErrorCode ierr
 13:       PetscInt indices(5),n
 14:       PetscInt five
 15:       PetscMPIInt rank
 16:       PetscInt, pointer :: idx(:)
 17:       IS      is

 19:       five = 5
 20:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 21:       if (ierr .ne. 0) then
 22:         print*,'Unable to initialize PETSc'
 23:         stop
 24:       endif
 25:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)

 27: !  Create an index set with 5 entries. Each processor creates
 28: !  its own index set with its own list of integers.

 30:       indices(1) = rank + 1
 31:       indices(2) = rank + 2
 32:       indices(3) = rank + 3
 33:       indices(4) = rank + 4
 34:       indices(5) = rank + 5
 35:       call ISCreateGeneral(PETSC_COMM_SELF,five,indices,PETSC_COPY_VALUES,is,ierr);CHKERRA(ierr)

 37: !  Print the index set to stdout

 39:       call ISView(is,PETSC_VIEWER_STDOUT_SELF,ierr);CHKERRA(ierr)

 41: !  Get the number of indices in the set

 43:       call ISGetLocalSize(is,n,ierr);CHKERRA(ierr)

 45: !   Get the indices in the index set

 47:       call ISGetIndicesF90(is,idx,ierr);CHKERRA(ierr)

 49:       if (associated(idx)) then
 50:          write (*,*) 'Association check passed'
 51:       else
 52:          write (*,*) 'Association check failed'
 53:       endif

 55: !   Now any code that needs access to the list of integers
 56: !   has access to it here

 58:       write(6,50) idx
 59:  50   format(5I3)

 61:       write(6,100) rank,idx(1),idx(5)
 62:  100  format('[',i5,'] First index = ',i5,' fifth index = ',i5)

 64: !   Once we no longer need access to the indices they should
 65: !   returned to the system

 67:       call ISRestoreIndicesF90(is,idx,ierr);CHKERRA(ierr)

 69: !   All PETSc objects should be destroyed once they are
 70: !   no longer needed

 72:       call ISDestroy(is,ierr);CHKERRA(ierr)
 73:       call PetscFinalize(ierr)
 74:       end

 76: !/*TEST
 77: !
 78: !   test:
 79: !      filter: sort -b
 80: !      filter_output: sort -b
 81: !
 82: !TEST*/