Actual source code: ex2f.F

  1: !
  2: !      "$Id: ex2f.F,v 1.19 2001/01/15 21:44:33 bsmith Exp $";
  3: !
  4: !    Description: Creates an index set based on a stride. Views that index set
  5: !    and then destroys it.
  6: !
  7: !/*T
  8: !    Concepts: index sets^manipulating a stride index set;
  9: !    Concepts: index sets^accessing indices from Fortran
 10: !T*/
 11: !
 12: !  Include petscis.h so we can use PETSc IS objects.
 13: !
 14:       program main
 15:       implicit none
 16:  #include finclude/petsc.h
 17:  #include finclude/petscis.h

 19:       integer     i,n,ierr,index(1),first,step,val
 20:       IS          set
 21:       PetscOffset iss

 23: #define indices(ib)  index(iss + (ib))

 25:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 26:       n     = 10
 27:       first = 3
 28:       step  = 2

 30: !  Create stride index set, starting at 3 with a stride of 2
 31: !    Note each processor is generating its own index set
 32: !    (in this case they are all identical)

 34:       call ISCreateStride(PETSC_COMM_SELF,n,first,step,set,ierr)
 35:       call ISView(set,PETSC_VIEWER_STDOUT_SELF,ierr)

 37: !  Extract the indice values from the set. Demonstrates how a Fortran code can
 38: !  directly access the array storing a PETSc index set with ISGetIndices().
 39: !  The user declares an array (index(1)) and index variable (iss), which are
 40: !  then used together to allow the Fortran to directly manipulate the PETSc array

 42:       call ISGetIndices(set,index,iss,ierr)
 43:       write(6,20)
 44:       do 10 i=1,n
 45: !  Bug in IRIX64 f90 compiler - write cannot handle integer(integer*8) correctly
 46:          val = indices(i)
 47:          write(6,30) val
 48:  10   continue
 49:  20   format('Printing indices directly')
 50:  30   format(i3)
 51:       call ISRestoreIndices(set,index,iss,ierr)

 53: !  Determine information on stride

 55:       call ISStrideGetInfo(set,first,step,ierr)
 56:       if (first .ne. 3 .or. step .ne. 2) then
 57:         print*,'Stride info not correct!'
 58:       endif

 60:       call ISDestroy(set,ierr)
 61:       call PetscFinalize(ierr)
 62:       end