Actual source code: ex1f.F90
1: !
2: !
3: ! Formatted test for IS general routines
4: !
5: #include <petsc/finclude/petscis.h>
6: program main
7: use petscis
8: implicit none
10: PetscErrorCode ierr
11: PetscInt i, n, indices(1004)
12: PetscInt, pointer :: ii(:)
13: PetscMPIInt size, rank
14: IS is, newis
15: PetscBool flag, compute, permanent
17: PetscCallA(PetscInitialize(ierr))
18: PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr))
19: PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD, size, ierr))
21: ! Test IS of size 0
22: n = 0
23: PetscCallA(ISCreateGeneral(PETSC_COMM_SELF, n, indices, PETSC_COPY_VALUES, is, ierr))
24: PetscCallA(ISGetLocalSize(is, n, ierr))
25: PetscCheckA(n == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, 'Error getting size of zero IS')
26: PetscCallA(ISDestroy(is, ierr))
28: ! Create large IS and test ISGetIndices(,ierr))
29: ! fortran indices start from 1 - but IS indices start from 0
30: n = 1000 + rank
31: do i = 1, n ! MarDiehl: will fail for rank > 4
32: indices(i) = rank + i - 1
33: end do
34: PetscCallA(ISCreateGeneral(PETSC_COMM_SELF, n, indices, PETSC_COPY_VALUES, is, ierr))
35: PetscCallA(ISGetIndices(is, ii, ierr))
36: do i = 1, n
37: PetscCheckA(ii(i) == indices(i), PETSC_COMM_SELF, PETSC_ERR_PLIB, 'Error getting indices')
38: end do
39: PetscCallA(ISRestoreIndices(is, ii, ierr))
41: ! Check identity and permutation
43: compute = PETSC_TRUE
44: permanent = PETSC_FALSE
45: PetscCallA(ISPermutation(is, flag, ierr))
46: PetscCheckA(.not. flag, PETSC_COMM_SELF, PETSC_ERR_PLIB, 'Error checking permutation')
47: PetscCallA(ISGetInfo(is, IS_PERMUTATION, IS_LOCAL, compute, flag, ierr))
48: PetscCallA(ISIdentity(is, flag, ierr))
49: PetscCheckA((rank /= 0) .or. flag, PETSC_COMM_SELF, PETSC_ERR_PLIB, 'Error checking identity')
50: PetscCheckA((rank == 0) .or. (.not. flag), PETSC_COMM_SELF, PETSC_ERR_PLIB, 'Error checking identity')
51: PetscCallA(ISSetInfo(is, IS_PERMUTATION, IS_LOCAL, permanent, PETSC_TRUE, ierr))
52: PetscCallA(ISSetInfo(is, IS_IDENTITY, IS_LOCAL, permanent, PETSC_TRUE, ierr))
53: PetscCallA(ISGetInfo(is, IS_PERMUTATION, IS_LOCAL, compute, flag, ierr))
54: PetscCheckA(flag, PETSC_COMM_SELF, PETSC_ERR_PLIB, 'Error checking permutation second time')
55: PetscCallA(ISGetInfo(is, IS_IDENTITY, IS_LOCAL, compute, flag, ierr))
56: PetscCheckA(flag, PETSC_COMM_SELF, PETSC_ERR_PLIB, 'Error checking identity second time')
57: PetscCallA(ISClearInfoCache(is, PETSC_TRUE, ierr))
59: ! Check equality of index sets
61: PetscCallA(ISEqual(is, is, flag, ierr))
62: PetscCheckA(flag, PETSC_COMM_SELF, PETSC_ERR_PLIB, 'Error checking equal')
64: ! Sorting
66: PetscCallA(ISSort(is, ierr))
67: PetscCallA(ISSorted(is, flag, ierr))
68: PetscCheckA(flag, PETSC_COMM_SELF, PETSC_ERR_PLIB, 'Error checking sorted')
70: ! Thinks it is a different type?
72: PetscCallA(PetscObjectTypeCompare(is, ISSTRIDE, flag, ierr))
73: PetscCheckA(.not. flag, PETSC_COMM_SELF, PETSC_ERR_PLIB, 'Error checking stride')
74: PetscCallA(PetscObjectTypeCompare(is, ISBLOCK, flag, ierr))
75: PetscCheckA(.not. flag, PETSC_COMM_SELF, PETSC_ERR_PLIB, 'Error checking block')
77: PetscCallA(ISDestroy(is, ierr))
79: ! Inverting permutation
81: do i = 1, n
82: indices(i) = n - i
83: end do
85: PetscCallA(ISCreateGeneral(PETSC_COMM_SELF, n, indices, PETSC_COPY_VALUES, is, ierr))
86: PetscCallA(ISSetPermutation(is, ierr))
87: PetscCallA(ISInvertPermutation(is, PETSC_DECIDE, newis, ierr))
88: PetscCallA(ISGetIndices(newis, ii, ierr))
89: do i = 1, n
90: PetscCheckA(ii(i) == n - i, PETSC_COMM_SELF, PETSC_ERR_PLIB, 'Error getting permutation indices')
91: end do
92: PetscCallA(ISRestoreIndices(newis, ii, ierr))
93: PetscCallA(ISDestroy(newis, ierr))
94: PetscCallA(ISDestroy(is, ierr))
95: PetscCallA(PetscFinalize(ierr))
96: end
98: !/*TEST
99: !
100: ! test:
101: ! nsize: {{1 2 3 4 5}}
102: ! output_file: output/empty.out
103: !
104: !TEST*/