xref: /petsc/src/sys/tests/ex41.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1c4762a1bSJed Brown static char help[] = "Test PETSc integer hash set.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petsc/private/hashseti.h>
4c4762a1bSJed Brown #include <petscsys.h>
5c4762a1bSJed Brown 
62c71b3e2SJacob Faibussowitsch #define PetscTestCheck(expr) PetscCheck(expr,PETSC_COMM_SELF,PETSC_ERR_LIB,"Assertion: `%s' failed.",PetscStringize(expr))
7c4762a1bSJed Brown 
8c4762a1bSJed Brown int main(int argc,char **argv)
9c4762a1bSJed Brown {
10c4762a1bSJed Brown   PetscHSetI     ht = NULL, hd;
11c4762a1bSJed Brown   PetscInt       n, off, array[4],na,nb,i,*marray,size;
12c4762a1bSJed Brown   PetscBool      has, flag;
13c4762a1bSJed Brown   PetscErrorCode ierr;
14c4762a1bSJed Brown 
15c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
16c4762a1bSJed Brown 
17*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetICreate(&ht));
182c71b3e2SJacob Faibussowitsch   PetscTestCheck(ht != NULL);
19*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIGetSize(ht,&n));
202c71b3e2SJacob Faibussowitsch   PetscTestCheck(n == 0);
21c4762a1bSJed Brown 
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIResize(ht,0));
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIGetSize(ht,&n));
242c71b3e2SJacob Faibussowitsch   PetscTestCheck(n == 0);
25c4762a1bSJed Brown 
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIHas(ht,42,&has));
272c71b3e2SJacob Faibussowitsch   PetscTestCheck(has == PETSC_FALSE);
28c4762a1bSJed Brown 
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIAdd(ht,42));
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIGetSize(ht,&n));
312c71b3e2SJacob Faibussowitsch   PetscTestCheck(n == 1);
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIHas(ht,42,&has));
332c71b3e2SJacob Faibussowitsch   PetscTestCheck(has == PETSC_TRUE);
34c4762a1bSJed Brown 
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIDel(ht,42));
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIGetSize(ht,&n));
372c71b3e2SJacob Faibussowitsch   PetscTestCheck(n == 0);
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIHas(ht,42,&has));
392c71b3e2SJacob Faibussowitsch   PetscTestCheck(has == PETSC_FALSE);
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIDel(ht,42));
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIDel(ht,24));
42c4762a1bSJed Brown 
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIQueryAdd(ht,123,&flag));
442c71b3e2SJacob Faibussowitsch   PetscTestCheck(flag == PETSC_TRUE);
45*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIQueryAdd(ht,123,&flag));
462c71b3e2SJacob Faibussowitsch   PetscTestCheck(flag == PETSC_FALSE);
47*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIQueryDel(ht,123,&flag));
482c71b3e2SJacob Faibussowitsch   PetscTestCheck(flag == PETSC_TRUE);
49*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIQueryDel(ht,123,&flag));
502c71b3e2SJacob Faibussowitsch   PetscTestCheck(flag == PETSC_FALSE);
51c4762a1bSJed Brown 
52*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIResize(ht,13));
53*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIGetSize(ht,&n));
542c71b3e2SJacob Faibussowitsch   PetscTestCheck(n == 0);
55c4762a1bSJed Brown 
56*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIClear(ht));
57*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIGetSize(ht,&n));
582c71b3e2SJacob Faibussowitsch   PetscTestCheck(n == 0);
59c4762a1bSJed Brown 
60*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIAdd(ht,42));
61*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIAdd(ht,13));
62*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIGetSize(ht,&n));
632c71b3e2SJacob Faibussowitsch   PetscTestCheck(n == 2);
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   off = 0;
66*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIGetElems(ht,&off,array));
67*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSortInt(off,array));
682c71b3e2SJacob Faibussowitsch   PetscTestCheck(off == 2);
692c71b3e2SJacob Faibussowitsch   PetscTestCheck(array[0] == 13);
702c71b3e2SJacob Faibussowitsch   PetscTestCheck(array[1] == 42);
71*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIGetElems(ht,&off,array));
72*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSortInt(2,array+2));
732c71b3e2SJacob Faibussowitsch   PetscTestCheck(off == 4);
742c71b3e2SJacob Faibussowitsch   PetscTestCheck(array[0] == 13);
752c71b3e2SJacob Faibussowitsch   PetscTestCheck(array[1] == 42);
762c71b3e2SJacob Faibussowitsch   PetscTestCheck(array[0] == 13);
772c71b3e2SJacob Faibussowitsch   PetscTestCheck(array[1] == 42);
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   off = 0;
80*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIDuplicate(ht,&hd));
81*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIGetElems(hd,&off,array));
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSortInt(off,array));
832c71b3e2SJacob Faibussowitsch   PetscTestCheck(off == 2);
842c71b3e2SJacob Faibussowitsch   PetscTestCheck(array[0] == 13);
852c71b3e2SJacob Faibussowitsch   PetscTestCheck(array[1] == 42);
86*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIDestroy(&hd));
87c4762a1bSJed Brown 
88*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIAdd(ht,0));
89*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIGetSize(ht,&n));
902c71b3e2SJacob Faibussowitsch   PetscTestCheck(n != 0);
91*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIReset(ht));
92*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIGetSize(ht,&n));
932c71b3e2SJacob Faibussowitsch   PetscTestCheck(n == 0);
94*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIReset(ht));
95*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIGetSize(ht,&n));
962c71b3e2SJacob Faibussowitsch   PetscTestCheck(n == 0);
97*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIAdd(ht,0));
98*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIGetSize(ht,&n));
992c71b3e2SJacob Faibussowitsch   PetscTestCheck(n != 0);
100c4762a1bSJed Brown 
101*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIDestroy(&ht));
1022c71b3e2SJacob Faibussowitsch   PetscTestCheck(ht == NULL);
103c4762a1bSJed Brown 
104*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetICreate(&ht));
105*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIReset(ht));
106*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIGetSize(ht,&n));
1072c71b3e2SJacob Faibussowitsch   PetscTestCheck(n == 0);
108*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIDestroy(&ht));
109c4762a1bSJed Brown 
110*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetICreate(&ht));
111*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetICreate(&hd));
112c4762a1bSJed Brown   n = 10;
113*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIResize(ht,n));
114*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIResize(hd,n));
115*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIGetCapacity(ht,&na));
116*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIGetCapacity(hd,&nb));
1172c71b3e2SJacob Faibussowitsch   PetscTestCheck(na>=n);
1182c71b3e2SJacob Faibussowitsch   PetscTestCheck(nb>=n);
119c4762a1bSJed Brown   for (i=0; i<n; i++) {
120*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscHSetIAdd(ht,i+1));
121*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscHSetIAdd(hd,i+1+n));
122c4762a1bSJed Brown   }
123*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIGetCapacity(ht,&nb));
1242c71b3e2SJacob Faibussowitsch   PetscTestCheck(nb>=na);
125c4762a1bSJed Brown   /* Merge ht and hd, and the result is in ht */
126*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIUpdate(ht,hd));
127*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIDestroy(&hd));
128*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIGetSize(ht,&size));
1292c71b3e2SJacob Faibussowitsch   PetscTestCheck(size==(2*n));
130*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n*2,&marray));
131c4762a1bSJed Brown   off = 0;
132*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIGetElems(ht,&off,marray));
133*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscHSetIDestroy(&ht));
1342c71b3e2SJacob Faibussowitsch   PetscTestCheck(off==(2*n));
135*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSortInt(off,marray));
136c4762a1bSJed Brown   for (i=0; i<n; i++) {
1372c71b3e2SJacob Faibussowitsch     PetscTestCheck(marray[i]==(i+1));
1382c71b3e2SJacob Faibussowitsch     PetscTestCheck(marray[n+i]==(i+1+n));
139c4762a1bSJed Brown   }
140*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(marray));
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   ierr = PetscFinalize();
143c4762a1bSJed Brown   return ierr;
144c4762a1bSJed Brown }
145c4762a1bSJed Brown 
146c4762a1bSJed Brown /*TEST
147c4762a1bSJed Brown 
148c4762a1bSJed Brown    test:
149c4762a1bSJed Brown 
150c4762a1bSJed Brown TEST*/
151