Actual source code: ex1f90.F90

petsc-3.10.5 2019-03-28
Report Typos and Errors
  1: !
  2: !  Description: Creates an index set based on a set of integers. Views that index set
  3: !  and then destroys it.
  4: !
  5: !/*T
  6: !    Concepts: index sets^manipulating a general index set;
  7: !    Concepts: Fortran90^accessing indices of index set;
  8: !T*/
  9: !

 11:       program main
 12:  #include <petsc/finclude/petscis.h>
 13:       use petscis
 14:       implicit none

 16:       PetscErrorCode ierr
 17:       PetscInt indices(5),n
 18:       PetscInt five
 19:       PetscMPIInt rank
 20:       PetscInt, pointer :: idx(:)
 21:       IS      is

 23:       five = 5
 24:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 25:       if (ierr .ne. 0) then
 26:         print*,'Unable to initialize PETSc'
 27:         stop
 28:       endif
 29:       call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)

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

 34:       indices(1) = rank + 1
 35:       indices(2) = rank + 2
 36:       indices(3) = rank + 3
 37:       indices(4) = rank + 4
 38:       indices(5) = rank + 5
 39:       call ISCreateGeneral(PETSC_COMM_SELF,five,indices,PETSC_COPY_VALUES,is,ierr);CHKERRA(ierr)

 41: !  Print the index set to stdout

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

 45: !  Get the number of indices in the set

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

 49: !   Get the indices in the index set

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

 53:       if (associated(idx)) then
 54:          write (*,*) 'Association check passed'
 55:       else
 56:          write (*,*) 'Association check failed'
 57:       endif

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

 62:       write(6,50) idx
 63:  50   format(5I3)

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

 68: !   Once we no longer need access to the indices they should
 69: !   returned to the system

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

 73: !   All PETSc objects should be destroyed once they are
 74: !   no longer needed

 76:       call ISDestroy(is,ierr);CHKERRA(ierr)
 77:       call PetscFinalize(ierr)
 78:       end

 80: !/*TEST
 81: !
 82: !   test:
 83: !      filter: sort -b
 84: !      filter_output: sort -b
 85: !
 86: !TEST*/