xref: /petsc/src/vec/is/is/tests/ex1.c (revision 6ffe77eaecce1557e50d00ca5347a7f48e598865)
1 /*
2        Formatted test for ISGeneral routines.
3 */
4 
5 static char help[] = "Tests IS general routines.\n\n";
6 
7 #include <petscis.h>
8 #include <petscviewer.h>
9 
10 int main(int argc,char **argv)
11 {
12   PetscMPIInt    rank,size;
13   PetscInt       i,n,*indices;
14   const PetscInt *ii;
15   IS             is,newis;
16   PetscBool      flg;
17   PetscBool      permanent = PETSC_FALSE;
18   PetscBool      compute = PETSC_TRUE;
19 
20   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
21   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
22   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
23 
24   /*
25      Test IS of size 0
26   */
27   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,0,&n,PETSC_COPY_VALUES,&is));
28   PetscCall(ISGetSize(is,&n));
29   PetscCheck(n == 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ISGetSize");
30   PetscCall(ISDestroy(&is));
31 
32   /*
33      Create large IS and test ISGetIndices()
34   */
35   n    = 10000 + rank;
36   PetscCall(PetscMalloc1(n,&indices));
37   for (i=0; i<n; i++) indices[i] = rank + i;
38   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,n,indices,PETSC_COPY_VALUES,&is));
39   PetscCall(ISGetIndices(is,&ii));
40   for (i=0; i<n; i++) {
41     PetscCheck(ii[i] == indices[i],PETSC_COMM_SELF,PETSC_ERR_PLIB,"ISGetIndices");
42   }
43   PetscCall(ISRestoreIndices(is,&ii));
44 
45   /*
46      Check identity and permutation
47   */
48   /* ISPermutation doesn't check if not set */
49   PetscCall(ISPermutation(is,&flg));
50   PetscCheck(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ISPermutation");
51   PetscCall(ISGetInfo(is,IS_PERMUTATION,IS_LOCAL,compute,&flg));
52   PetscCheck(rank != 0 || flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ISGetInfo(IS_PERMUTATION,IS_LOCAL)");
53   PetscCheck(rank == 0 || !flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ISGetInfo(IS_PERMUTATION,IS_LOCAL)");
54   PetscCall(ISIdentity(is,&flg));
55   PetscCheck(rank != 0 || flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ISIdentity");
56   PetscCheck(rank == 0 || !flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ISIdentity");
57   PetscCall(ISGetInfo(is,IS_IDENTITY,IS_LOCAL,compute,&flg));
58   PetscCheck(rank != 0 || flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ISGetInfo(IS_IDENTITY,IS_LOCAL)");
59   PetscCheck(rank == 0 || !flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ISGetInfo(IS_IDENTITY,IS_LOCAL)");
60   /* we can override the computed values with ISSetInfo() */
61   PetscCall(ISSetInfo(is,IS_PERMUTATION,IS_LOCAL,permanent,PETSC_TRUE));
62   PetscCall(ISSetInfo(is,IS_IDENTITY,IS_LOCAL,permanent,PETSC_TRUE));
63   PetscCall(ISGetInfo(is,IS_PERMUTATION,IS_LOCAL,compute,&flg));
64   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ISGetInfo(IS_PERMUTATION,IS_LOCAL)");
65   PetscCall(ISGetInfo(is,IS_IDENTITY,IS_LOCAL,compute,&flg));
66   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ISGetInfo(IS_IDENTITY,IS_LOCAL)");
67 
68   PetscCall(ISClearInfoCache(is,PETSC_TRUE));
69 
70   /*
71      Check equality of index sets
72   */
73   PetscCall(ISEqual(is,is,&flg));
74   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ISEqual");
75 
76   /*
77      Sorting
78   */
79   PetscCall(ISSort(is));
80   PetscCall(ISSorted(is,&flg));
81   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ISSort");
82   PetscCall(ISGetInfo(is,IS_SORTED,IS_LOCAL,compute,&flg));
83   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ISGetInfo(IS_SORTED,IS_LOCAL)");
84   PetscCall(ISSorted(is,&flg));
85   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ISSort");
86   PetscCall(ISGetInfo(is,IS_SORTED,IS_LOCAL,compute,&flg));
87   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ISGetInfo(IS_SORTED,IS_LOCAL)");
88 
89   /*
90      Thinks it is a different type?
91   */
92   PetscCall(PetscObjectTypeCompare((PetscObject)is,ISSTRIDE,&flg));
93   PetscCheck(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ISStride");
94   PetscCall(PetscObjectTypeCompare((PetscObject)is,ISBLOCK,&flg));
95   PetscCheck(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ISBlock");
96 
97   PetscCall(ISDestroy(&is));
98 
99   /*
100      Inverting permutation
101   */
102   for (i=0; i<n; i++) indices[i] = n - i - 1;
103   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,n,indices,PETSC_COPY_VALUES,&is));
104   PetscCall(PetscFree(indices));
105   PetscCall(ISSetPermutation(is));
106   PetscCall(ISInvertPermutation(is,PETSC_DECIDE,&newis));
107   PetscCall(ISGetIndices(newis,&ii));
108   for (i=0; i<n; i++) {
109     PetscCheck(ii[i] == n - i - 1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ISInvertPermutation");
110   }
111   PetscCall(ISRestoreIndices(newis,&ii));
112   PetscCall(ISDestroy(&newis));
113   PetscCall(ISDestroy(&is));
114   PetscCall(PetscFinalize());
115   return 0;
116 }
117 
118 /*TEST
119 
120    test:
121       nsize: {{1 2 3 4 5}}
122 
123 TEST*/
124