Actual source code: ex2f.F

petsc-3.10.5 2019-03-28
Report Typos and Errors
  1: !
  2: !     Description: Creates an index set based on a stride. Views that
  3: !     index set and then destroys it.
  4: !
  5: !/*T
  6: !     Concepts: index sets^manipulating a stride index set;
  7: !     Concepts: index sets^accessing indices from Fortran
  8: !T*/
  9: !
 10: !     Include petscis.h so we can use PETSc IS objects.
 11: !
 12:       program main
 13:  #include <petsc/finclude/petscis.h>
 14:       use petscis
 15:       implicit none

 17:       PetscErrorCode ierr
 18:       PetscInt    i,n,index(1),first,step,val
 19:       IS          set
 20:       PetscOffset iss

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

 24:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 25:       if (ierr .ne. 0) then
 26:          print*,'Unable to initialize PETSc'
 27:          stop
 28:       endif
 29:       n     = 10
 30:       first = 3
 31:       step  = 2

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

 37:       call ISCreateStride(PETSC_COMM_SELF,n,first,step,set,ierr)
 38:       call ISView(set,PETSC_VIEWER_STDOUT_SELF,ierr)

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

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

 58: !     Determine information on stride

 60:       call ISStrideGetInfo(set,first,step,ierr)
 61:       if (first .ne. 3 .or. step .ne. 2) then
 62:         print*,'Stride info not correct!'
 63:       endif

 65:       call ISDestroy(set,ierr)
 66:       call PetscFinalize(ierr)
 67:       end


 70: !/*TEST
 71: !
 72: !   test:
 73: !
 74: !TEST*/