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: !
  6: #include <petsc/finclude/petscis.h>
  7: program main
  8:   use petscis
  9:   implicit none

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

 17:   PetscCallA(PetscInitialize(ierr))
 18:   PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr))

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

 23:   do i = 1, 5_PETSC_INT_KIND
 24:     indices(i) = rank + i
 25:   end do
 26:   PetscCallA(ISCreateGeneral(PETSC_COMM_SELF, 5_PETSC_INT_KIND, indices, PETSC_COPY_VALUES, is, ierr))

 28: !  Print the index set to stdout

 30:   PetscCallA(ISView(is, PETSC_VIEWER_STDOUT_SELF, ierr))

 32: !  Get the number of indices in the set

 34:   PetscCallA(ISGetLocalSize(is, n, ierr))

 36: !   Get the indices in the index set

 38:   PetscCallA(ISGetIndices(is, idx, ierr))

 40:   if (associated(idx)) then
 41:     write (*, *) 'Association check passed'
 42:   else
 43:     write (*, *) 'Association check failed'
 44:   end if

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

 49:   write (6, 50) idx
 50: 50 format(5I3)

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

 55: !   Once we no longer need access to the indices they should
 56: !   returned to the system

 58:   PetscCallA(ISRestoreIndices(is, idx, ierr))

 60: !   All PETSc objects should be destroyed once they are
 61: !   no longer needed

 63:   PetscCallA(ISDestroy(is, ierr))
 64:   PetscCallA(PetscFinalize(ierr))
 65: end

 67: !/*TEST
 68: !
 69: !   test:
 70: !      filter: sort -b
 71: !      filter_output: sort -b
 72: !
 73: !TEST*/