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: #include <petsc/finclude/petscis.h>
9: program main
10: use petscis
11: implicit none
13: PetscErrorCode ierr
14: PetscInt i, first, step
15: PetscInt, parameter :: n = 10
16: IS set
17: PetscInt, pointer :: index(:)
19: PetscCallA(PetscInitialize(ierr))
20: first = 3
21: step = 2
23: ! Create stride index set, starting at 3 with a stride of 2 Note
24: ! each processor is generating its own index set (in this case they
25: ! are all identical)
27: PetscCallA(ISCreateStride(PETSC_COMM_SELF, n, first, step, set, ierr))
28: PetscCallA(ISView(set, PETSC_VIEWER_STDOUT_SELF, ierr))
30: ! Extract the indices values from the set. Demonstrates how a Fortran
31: ! code can directly access the array storing a PETSc index set with
32: ! ISGetIndices().
34: PetscCallA(ISGetIndices(set, index, ierr))
35: write (6, 20)
36: do i = 1, n
37: write (6, 30) index(i)
38: end do
39: 20 format('Printing indices directly')
40: 30 format(i3)
41: PetscCallA(ISRestoreIndices(set, index, ierr))
43: ! Determine information on stride
45: PetscCallA(ISStrideGetInfo(set, first, step, ierr))
46: if (first /= 3 .or. step /= 2) print *, 'Stride info not correct!'
48: PetscCallA(ISDestroy(set, ierr))
49: PetscCallA(PetscFinalize(ierr))
50: end
52: !/*TEST
53: !
54: ! test:
55: !
56: !TEST*/