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