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