Actual source code: ex3f90.F

  1: !
  2: !    Description:  Creates an index set based on blocks of integers. Views that index set
  3: !    and then destroys it.
  4: !
  5: !/*T
  6: !    Concepts: index sets^manipulating a block index set;
  7: !    Concepts: Fortran90^accessing indices in index set;
  8: !
  9: !T*/
 10: !
 11: !  The following include statements are required for Fortran programs
 12: !  that use PETSc index sets:
 13: !     petsc.h  - base PETSc routines
 14: !     petscis.h     - index sets (IS objects)
 15: !     petscis.h90   - to allow access to Fortran90 features of index sets
 16: !
 17:       program main
 18:       implicit none

 20:  #include finclude/petsc.h
 21:  #include finclude/petscis.h
 22: #include "finclude/petscis.h90"

 24:       integer n,ierr,bs,issize,inputindices(4)
 25:       integer, pointer :: indices(:)
 26:       IS       set
 27:       PetscTruth isablock;

 29:       n               = 4
 30:       bs              = 3
 31:       inputindices(1) = 0
 32:       inputindices(2) = 3
 33:       inputindices(3) = 9
 34:       inputindices(4) = 12
 35: 
 36:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 37: 
 38: !
 39: !    Create a block index set. The index set has 4 blocks each of size 3.
 40: !    The indices are {0,1,2,3,4,5,9,10,11,12,13,14}
 41: !    Note each processor is generating its own index set
 42: !    (in this case they are all identical)
 43: !
 44:       call ISCreateBlock(PETSC_COMM_SELF,bs,n,inputindices,set,ierr)
 45:       call ISView(set,PETSC_VIEWER_STDOUT_SELF,ierr)

 47: !
 48: !    Extract indices from set.
 49: !
 50:       call ISGetLocalSize(set,issize,ierr)
 51:       call ISGetIndicesF90(set,indices,ierr)
 52:       write(6,100)indices
 53:  100  format(12I3)
 54:       call ISRestoreIndicesF90(set,indices,ierr)

 56: !
 57: !    Extract the block indices. This returns one index per block.
 58: !
 59:       call ISBlockGetIndicesF90(set,indices,ierr)
 60:       write(6,200)indices
 61:  200  format(4I3)
 62:       call ISBlockRestoreIndicesF90(set,indices,ierr)

 64: !
 65: !    Check if this is really a block index set
 66: !
 67:       call ISBlock(set,isablock,ierr)
 68:       if (isablock .ne. PETSC_TRUE) then
 69:         write(6,*) 'Index set is not blocked!'
 70:       endif

 72: !
 73: !    Determine the block size of the index set
 74: !
 75:       call ISBlockGetBlockSize(set,bs,ierr)
 76:       if (bs .ne. 3) then
 77:         write(6,*) 'Blocksize != 3'
 78:       endif

 80: !
 81: !    Get the number of blocks
 82: !
 83:       call ISBlockGetSize(set,n,ierr)
 84:       if (n .ne. 4) then
 85:         write(6,*) 'Number of blocks != 4'
 86:       endif

 88:       call ISDestroy(set,ierr)
 89:       call PetscFinalize(ierr)
 90:       end