Actual source code: ex2f.F90
1: !
2: ! Formatted Test for IS stride routines
3: !
4: #include <petsc/finclude/petscis.h>
5: program main
6: use petscis
7: implicit none
9: PetscErrorCode ierr
10: PetscInt i, n, start
11: PetscInt stride, ssize, first
12: IS is
13: PetscBool flag
14: PetscInt, pointer :: ii(:)
16: PetscCallA(PetscInitialize(ierr))
18: ! Test IS of size 0
19: ssize = 0
20: stride = 0
21: first = 2
22: PetscCallA(ISCreateStride(PETSC_COMM_SELF, ssize, stride, first, is, ierr))
23: PetscCallA(ISGetLocalSize(is, n, ierr))
24: PetscCheckA(n == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, 'Wrong result from ISCreateStride')
26: PetscCallA(ISStrideGetInfo(is, start, stride, ierr))
27: PetscCheckA(start == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, 'Wrong result from ISStrideGetInfo')
28: PetscCheckA(stride == 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, 'Wrong result from ISStrideGetInfo')
30: PetscCallA(PetscObjectTypeCompare(is, ISSTRIDE, flag, ierr))
31: PetscCheckA(flag, PETSC_COMM_SELF, PETSC_ERR_PLIB, 'Wrong result from PetscObjectTypeCompare')
32: PetscCallA(ISGetIndices(is, ii, ierr))
33: PetscCallA(ISRestoreIndices(is, ii, ierr))
34: PetscCallA(ISDestroy(is, ierr))
36: ! Test ISGetIndices()
37: ssize = 10000
38: stride = -8
39: first = 3
40: PetscCallA(ISCreateStride(PETSC_COMM_SELF, ssize, stride, first, is, ierr))
41: PetscCallA(ISGetLocalSize(is, n, ierr))
42: PetscCallA(ISGetIndices(is, ii, ierr))
43: do i = 1, n
44: PetscCheckA(ii(i) == -11 + 3*i, PETSC_COMM_SELF, PETSC_ERR_PLIB, 'Wrong result from ISGetIndices')
45: end do
46: PetscCallA(ISRestoreIndices(is, ii, ierr))
47: PetscCallA(ISDestroy(is, ierr))
49: PetscCallA(PetscFinalize(ierr))
50: end
52: !/*TEST
53: !
54: ! test:
55: ! output_file: output/empty.out
56: !
57: !TEST*/