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:   PetscCallA(PetscInitialize(ierr))
 21:   PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr))

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

 26:   indices(1) = rank + 1
 27:   indices(2) = rank + 2
 28:   indices(3) = rank + 3
 29:   indices(4) = rank + 4
 30:   indices(5) = rank + 5
 31:   PetscCallA(ISCreateGeneral(PETSC_COMM_SELF, five, indices, PETSC_COPY_VALUES, is, ierr))

 33: !  Print the index set to stdout

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

 37: !  Get the number of indices in the set

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

 41: !   Get the indices in the index set

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

 45:   if (associated(idx)) then
 46:     write (*, *) 'Association check passed'
 47:   else
 48:     write (*, *) 'Association check failed'
 49:   end if

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

 54:   write (6, 50) idx
 55: 50 format(5I3)

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

 60: !   Once we no longer need access to the indices they should
 61: !   returned to the system

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

 65: !   All PETSc objects should be destroyed once they are
 66: !   no longer needed

 68:   PetscCallA(ISDestroy(is, ierr))
 69:   PetscCallA(PetscFinalize(ierr))
 70: end

 72: !/*TEST
 73: !
 74: !   test:
 75: !      filter: sort -b
 76: !      filter_output: sort -b
 77: !
 78: !TEST*/