Actual source code: ex2f.F90
1: !
2: ! Description: Creates an index set based on a stride. Views that
3: ! index set and then destroys it.
4: !
5: !
6: ! Include petscis.h so we can use PETSc IS objects.
7: !
8: program main
9: #include <petsc/finclude/petscis.h>
10: use petscis
11: implicit none
13: PetscErrorCode ierr
14: PetscInt i,n,index(1),first,step,val
15: IS set
16: PetscOffset iss
18: #define indices(ib) index(iss + (ib))
20: PetscCallA(PetscInitialize(ierr))
21: n = 10
22: first = 3
23: step = 2
25: ! Create stride index set, starting at 3 with a stride of 2 Note
26: ! each processor is generating its own index set (in this case they
27: ! are all identical)
29: PetscCallA(ISCreateStride(PETSC_COMM_SELF,n,first,step,set,ierr))
30: PetscCallA(ISView(set,PETSC_VIEWER_STDOUT_SELF,ierr))
32: ! Extract the indice values from the set. Demonstrates how a Fortran
33: ! code can directly access the array storing a PETSc index set with
34: ! ISGetIndices(). The user declares an array (index(1)) and index
35: ! variable (iss), which are then used together to allow the Fortran
36: ! to directly manipulate the PETSc array
38: PetscCallA(ISGetIndices(set,index,iss,ierr))
39: write(6,20)
40: ! Bug in IRIX64 f90 compiler - write cannot handle
41: ! integer(integer*8) correctly
42: do 10 i=1,n
43: val = indices(i)
44: write(6,30) val
45: 10 continue
46: 20 format('Printing indices directly')
47: 30 format(i3)
48: PetscCallA(ISRestoreIndices(set,index,iss,ierr))
50: ! Determine information on stride
52: PetscCallA(ISStrideGetInfo(set,first,step,ierr))
53: if (first .ne. 3 .or. step .ne. 2) then
54: print*,'Stride info not correct!'
55: endif
57: PetscCallA(ISDestroy(set,ierr))
58: PetscCallA(PetscFinalize(ierr))
59: end
61: !/*TEST
62: !
63: ! test:
64: !
65: !TEST*/