Actual source code: ex1f.F90
1: !
2: !
3: ! Formatted test for IS general routines
4: !
5: program main
6: #include <petsc/finclude/petscis.h>
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
23: n = 0
24: PetscCallA(ISCreateGeneral(PETSC_COMM_SELF,n,indices,PETSC_COPY_VALUES,is,ierr))
25: PetscCallA(ISGetLocalSize(is,n,ierr))
26: PetscCheckA(n .eq. 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,'Error getting size of zero IS')
27: PetscCallA(ISDestroy(is,ierr))
29: ! Create large IS and test ISGetIndices(,ierr))
30: ! fortran indices start from 1 - but IS indices start from 0
31: n = 1000 + rank
32: do 10, i=1,n
33: indices(i) = rank + i-1
34: 10 continue
35: PetscCallA(ISCreateGeneral(PETSC_COMM_SELF,n,indices,PETSC_COPY_VALUES,is,ierr))
36: PetscCallA(ISGetIndicesF90(is,ii,ierr))
37: do 20, i=1,n
38: PetscCheckA(ii(i) .eq. indices(i),PETSC_COMM_SELF,PETSC_ERR_PLIB,'Error getting indices')
39: 20 continue
40: PetscCallA(ISRestoreIndicesF90(is,ii,ierr))
42: ! Check identity and permutation
44: compute = PETSC_TRUE
45: permanent = PETSC_FALSE
46: PetscCallA(ISPermutation(is,flag,ierr))
47: PetscCheckA(.not. flag,PETSC_COMM_SELF,PETSC_ERR_PLIB,'Error checking permutation')
48: PetscCallA(ISGetInfo(is,IS_PERMUTATION,IS_LOCAL,compute,flag,ierr))
49: PetscCallA(ISIdentity(is,flag,ierr))
50: PetscCheckA((rank .ne. 0) .or. flag,PETSC_COMM_SELF,PETSC_ERR_PLIB,'Error checking identity')
51: PetscCheckA((rank .eq. 0) .or. (.not. flag),PETSC_COMM_SELF,PETSC_ERR_PLIB,'Error checking identity')
52: PetscCallA(ISSetInfo(is,IS_PERMUTATION,IS_LOCAL,permanent,PETSC_TRUE,ierr))
53: PetscCallA(ISSetInfo(is,IS_IDENTITY,IS_LOCAL,permanent,PETSC_TRUE,ierr))
54: PetscCallA(ISGetInfo(is,IS_PERMUTATION,IS_LOCAL,compute,flag,ierr))
55: PetscCheckA(flag,PETSC_COMM_SELF,PETSC_ERR_PLIB,'Error checking permutation second time')
56: PetscCallA(ISGetInfo(is,IS_IDENTITY,IS_LOCAL,compute,flag,ierr))
57: PetscCheckA(flag,PETSC_COMM_SELF,PETSC_ERR_PLIB,'Error checking identity second time')
58: PetscCallA(ISClearInfoCache(is,PETSC_TRUE,ierr))
60: ! Check equality of index sets
62: PetscCallA(ISEqual(is,is,flag,ierr))
63: PetscCheckA(flag,PETSC_COMM_SELF,PETSC_ERR_PLIB,'Error checking equal')
65: ! Sorting
67: PetscCallA(ISSort(is,ierr))
68: PetscCallA(ISSorted(is,flag,ierr))
69: PetscCheckA(flag,PETSC_COMM_SELF,PETSC_ERR_PLIB,'Error checking sorted')
71: ! Thinks it is a different type?
73: PetscCallA(PetscObjectTypeCompare(is,ISSTRIDE,flag,ierr))
74: PetscCheckA(.not. flag,PETSC_COMM_SELF,PETSC_ERR_PLIB,'Error checking stride')
75: PetscCallA(PetscObjectTypeCompare(is,ISBLOCK,flag,ierr))
76: PetscCheckA(.not. flag,PETSC_COMM_SELF,PETSC_ERR_PLIB,'Error checking block')
78: PetscCallA(ISDestroy(is,ierr))
80: ! Inverting permutation
82: do 30, i=1,n
83: indices(i) = n - i
84: 30 continue
86: PetscCallA(ISCreateGeneral(PETSC_COMM_SELF,n,indices,PETSC_COPY_VALUES,is,ierr))
87: PetscCallA(ISSetPermutation(is,ierr))
88: PetscCallA(ISInvertPermutation(is,PETSC_DECIDE,newis,ierr))
89: PetscCallA(ISGetIndicesF90(newis,ii,ierr))
90: do 40, i=1,n
91: PetscCheckA(ii(i) .eq. n - i,PETSC_COMM_SELF,PETSC_ERR_PLIB,'Error getting permutation indices')
92: 40 continue
93: PetscCallA(ISRestoreIndicesF90(newis,ii,ierr))
94: PetscCallA(ISDestroy(newis,ierr))
95: PetscCallA(ISDestroy(is,ierr))
96: PetscCallA(PetscFinalize(ierr))
97: end
99: !/*TEST
100: !
101: ! test:
102: ! nsize: {{1 2 3 4 5}}
103: ! output_file: output/ex1_1.out
104: !
105: !TEST*/