xref: /petsc/src/sys/tests/ex41.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Test PETSc integer hash set.\n\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown #include <petsc/private/hashseti.h>
4*c4762a1bSJed Brown #include <petscsys.h>
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown #define PetscAssert(expr) do {            \
7*c4762a1bSJed Brown if (PetscUnlikely(!(expr)))               \
8*c4762a1bSJed Brown   SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, \
9*c4762a1bSJed Brown            "Assertion: `%s' failed.",     \
10*c4762a1bSJed Brown            PetscStringize(expr));         \
11*c4762a1bSJed Brown } while(0)
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown 
14*c4762a1bSJed Brown int main(int argc,char **argv)
15*c4762a1bSJed Brown {
16*c4762a1bSJed Brown   PetscHSetI     ht = NULL, hd;
17*c4762a1bSJed Brown   PetscInt       n, off, array[4],na,nb,i,*marray,size;
18*c4762a1bSJed Brown   PetscBool      has, flag;
19*c4762a1bSJed Brown   PetscErrorCode ierr;
20*c4762a1bSJed Brown 
21*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
22*c4762a1bSJed Brown 
23*c4762a1bSJed Brown   ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
24*c4762a1bSJed Brown   PetscAssert(ht != NULL);
25*c4762a1bSJed Brown   ierr = PetscHSetIGetSize(ht,&n);CHKERRQ(ierr);
26*c4762a1bSJed Brown   PetscAssert(n == 0);
27*c4762a1bSJed Brown 
28*c4762a1bSJed Brown   ierr = PetscHSetIResize(ht,0);CHKERRQ(ierr);
29*c4762a1bSJed Brown   ierr = PetscHSetIGetSize(ht,&n);CHKERRQ(ierr);
30*c4762a1bSJed Brown   PetscAssert(n == 0);
31*c4762a1bSJed Brown 
32*c4762a1bSJed Brown   ierr = PetscHSetIHas(ht,42,&has);CHKERRQ(ierr);
33*c4762a1bSJed Brown   PetscAssert(has == PETSC_FALSE);
34*c4762a1bSJed Brown 
35*c4762a1bSJed Brown   ierr = PetscHSetIAdd(ht,42);CHKERRQ(ierr);
36*c4762a1bSJed Brown   ierr = PetscHSetIGetSize(ht,&n);CHKERRQ(ierr);
37*c4762a1bSJed Brown   PetscAssert(n == 1);
38*c4762a1bSJed Brown   ierr = PetscHSetIHas(ht,42,&has);CHKERRQ(ierr);
39*c4762a1bSJed Brown   PetscAssert(has == PETSC_TRUE);
40*c4762a1bSJed Brown 
41*c4762a1bSJed Brown   ierr = PetscHSetIDel(ht,42);CHKERRQ(ierr);
42*c4762a1bSJed Brown   ierr = PetscHSetIGetSize(ht,&n);CHKERRQ(ierr);
43*c4762a1bSJed Brown   PetscAssert(n == 0);
44*c4762a1bSJed Brown   ierr = PetscHSetIHas(ht,42,&has);CHKERRQ(ierr);
45*c4762a1bSJed Brown   PetscAssert(has == PETSC_FALSE);
46*c4762a1bSJed Brown   ierr = PetscHSetIDel(ht,42);CHKERRQ(ierr);
47*c4762a1bSJed Brown   ierr = PetscHSetIDel(ht,24);CHKERRQ(ierr);
48*c4762a1bSJed Brown 
49*c4762a1bSJed Brown   ierr = PetscHSetIQueryAdd(ht,123,&flag);CHKERRQ(ierr);
50*c4762a1bSJed Brown   PetscAssert(flag == PETSC_TRUE);
51*c4762a1bSJed Brown   ierr = PetscHSetIQueryAdd(ht,123,&flag);CHKERRQ(ierr);
52*c4762a1bSJed Brown   PetscAssert(flag == PETSC_FALSE);
53*c4762a1bSJed Brown   ierr = PetscHSetIQueryDel(ht,123,&flag);CHKERRQ(ierr);
54*c4762a1bSJed Brown   PetscAssert(flag == PETSC_TRUE);
55*c4762a1bSJed Brown   ierr = PetscHSetIQueryDel(ht,123,&flag);CHKERRQ(ierr);
56*c4762a1bSJed Brown   PetscAssert(flag == PETSC_FALSE);
57*c4762a1bSJed Brown 
58*c4762a1bSJed Brown   ierr = PetscHSetIResize(ht,13);CHKERRQ(ierr);
59*c4762a1bSJed Brown   ierr = PetscHSetIGetSize(ht,&n);CHKERRQ(ierr);
60*c4762a1bSJed Brown   PetscAssert(n == 0);
61*c4762a1bSJed Brown 
62*c4762a1bSJed Brown   ierr = PetscHSetIClear(ht);CHKERRQ(ierr);
63*c4762a1bSJed Brown   ierr = PetscHSetIGetSize(ht,&n);CHKERRQ(ierr);
64*c4762a1bSJed Brown   PetscAssert(n == 0);
65*c4762a1bSJed Brown 
66*c4762a1bSJed Brown   ierr = PetscHSetIAdd(ht,42);CHKERRQ(ierr);
67*c4762a1bSJed Brown   ierr = PetscHSetIAdd(ht,13);CHKERRQ(ierr);
68*c4762a1bSJed Brown   ierr = PetscHSetIGetSize(ht,&n);CHKERRQ(ierr);
69*c4762a1bSJed Brown   PetscAssert(n == 2);
70*c4762a1bSJed Brown 
71*c4762a1bSJed Brown   off = 0;
72*c4762a1bSJed Brown   ierr = PetscHSetIGetElems(ht,&off,array);CHKERRQ(ierr);
73*c4762a1bSJed Brown   ierr = PetscSortInt(off,array);CHKERRQ(ierr);
74*c4762a1bSJed Brown   PetscAssert(off == 2);
75*c4762a1bSJed Brown   PetscAssert(array[0] == 13);
76*c4762a1bSJed Brown   PetscAssert(array[1] == 42);
77*c4762a1bSJed Brown   ierr = PetscHSetIGetElems(ht,&off,array);CHKERRQ(ierr);
78*c4762a1bSJed Brown   ierr = PetscSortInt(2,array+2);CHKERRQ(ierr);
79*c4762a1bSJed Brown   PetscAssert(off == 4);
80*c4762a1bSJed Brown   PetscAssert(array[0] == 13);
81*c4762a1bSJed Brown   PetscAssert(array[1] == 42);
82*c4762a1bSJed Brown   PetscAssert(array[0] == 13);
83*c4762a1bSJed Brown   PetscAssert(array[1] == 42);
84*c4762a1bSJed Brown 
85*c4762a1bSJed Brown   off = 0;
86*c4762a1bSJed Brown   ierr = PetscHSetIDuplicate(ht,&hd);CHKERRQ(ierr);
87*c4762a1bSJed Brown   ierr = PetscHSetIGetElems(hd,&off,array);CHKERRQ(ierr);
88*c4762a1bSJed Brown   ierr = PetscSortInt(off,array);CHKERRQ(ierr);
89*c4762a1bSJed Brown   PetscAssert(off == 2);
90*c4762a1bSJed Brown   PetscAssert(array[0] == 13);
91*c4762a1bSJed Brown   PetscAssert(array[1] == 42);
92*c4762a1bSJed Brown   ierr = PetscHSetIDestroy(&hd);CHKERRQ(ierr);
93*c4762a1bSJed Brown 
94*c4762a1bSJed Brown   ierr = PetscHSetIAdd(ht,0);CHKERRQ(ierr);
95*c4762a1bSJed Brown   ierr = PetscHSetIGetSize(ht,&n);CHKERRQ(ierr);
96*c4762a1bSJed Brown   PetscAssert(n != 0);
97*c4762a1bSJed Brown   ierr = PetscHSetIReset(ht);CHKERRQ(ierr);
98*c4762a1bSJed Brown   ierr = PetscHSetIGetSize(ht,&n);CHKERRQ(ierr);
99*c4762a1bSJed Brown   PetscAssert(n == 0);
100*c4762a1bSJed Brown   ierr = PetscHSetIReset(ht);CHKERRQ(ierr);
101*c4762a1bSJed Brown   ierr = PetscHSetIGetSize(ht,&n);CHKERRQ(ierr);
102*c4762a1bSJed Brown   PetscAssert(n == 0);
103*c4762a1bSJed Brown   ierr = PetscHSetIAdd(ht,0);CHKERRQ(ierr);
104*c4762a1bSJed Brown   ierr = PetscHSetIGetSize(ht,&n);CHKERRQ(ierr);
105*c4762a1bSJed Brown   PetscAssert(n != 0);
106*c4762a1bSJed Brown 
107*c4762a1bSJed Brown   ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr);
108*c4762a1bSJed Brown   PetscAssert(ht == NULL);
109*c4762a1bSJed Brown 
110*c4762a1bSJed Brown   ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
111*c4762a1bSJed Brown   ierr = PetscHSetIReset(ht);CHKERRQ(ierr);
112*c4762a1bSJed Brown   ierr = PetscHSetIGetSize(ht,&n);CHKERRQ(ierr);
113*c4762a1bSJed Brown   PetscAssert(n == 0);
114*c4762a1bSJed Brown   ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr);
115*c4762a1bSJed Brown 
116*c4762a1bSJed Brown   ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
117*c4762a1bSJed Brown   ierr = PetscHSetICreate(&hd);CHKERRQ(ierr);
118*c4762a1bSJed Brown   n = 10;
119*c4762a1bSJed Brown   ierr = PetscHSetIResize(ht,n);CHKERRQ(ierr);
120*c4762a1bSJed Brown   ierr = PetscHSetIResize(hd,n);CHKERRQ(ierr);
121*c4762a1bSJed Brown   ierr = PetscHSetIGetCapacity(ht,&na);CHKERRQ(ierr);
122*c4762a1bSJed Brown   ierr = PetscHSetIGetCapacity(hd,&nb);CHKERRQ(ierr);
123*c4762a1bSJed Brown   PetscAssert(na>=n);
124*c4762a1bSJed Brown   PetscAssert(nb>=n);
125*c4762a1bSJed Brown   for (i=0; i<n; i++) {
126*c4762a1bSJed Brown     ierr = PetscHSetIAdd(ht,i+1);CHKERRQ(ierr);
127*c4762a1bSJed Brown     ierr = PetscHSetIAdd(hd,i+1+n);CHKERRQ(ierr);
128*c4762a1bSJed Brown   }
129*c4762a1bSJed Brown   ierr = PetscHSetIGetCapacity(ht,&nb);CHKERRQ(ierr);
130*c4762a1bSJed Brown   PetscAssert(nb>=na);
131*c4762a1bSJed Brown   /* Merge ht and hd, and the result is in ht */
132*c4762a1bSJed Brown   ierr = PetscHSetIUpdate(ht,hd);CHKERRQ(ierr);
133*c4762a1bSJed Brown   ierr = PetscHSetIDestroy(&hd);CHKERRQ(ierr);
134*c4762a1bSJed Brown   ierr = PetscHSetIGetSize(ht,&size);CHKERRQ(ierr);
135*c4762a1bSJed Brown   PetscAssert(size==(2*n));CHKERRQ(ierr);
136*c4762a1bSJed Brown   ierr = PetscMalloc1(n*2,&marray);CHKERRQ(ierr);
137*c4762a1bSJed Brown   off = 0;
138*c4762a1bSJed Brown   ierr = PetscHSetIGetElems(ht,&off,marray);CHKERRQ(ierr);
139*c4762a1bSJed Brown   ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr);
140*c4762a1bSJed Brown   PetscAssert(off==(2*n));
141*c4762a1bSJed Brown   ierr = PetscSortInt(off,marray);CHKERRQ(ierr);
142*c4762a1bSJed Brown   for (i=0; i<n; i++) {
143*c4762a1bSJed Brown     PetscAssert(marray[i]==(i+1));
144*c4762a1bSJed Brown     PetscAssert(marray[n+i]==(i+1+n));
145*c4762a1bSJed Brown   }
146*c4762a1bSJed Brown   ierr = PetscFree(marray);CHKERRQ(ierr);
147*c4762a1bSJed Brown 
148*c4762a1bSJed Brown   ierr = PetscFinalize();
149*c4762a1bSJed Brown   return ierr;
150*c4762a1bSJed Brown }
151*c4762a1bSJed Brown 
152*c4762a1bSJed Brown 
153*c4762a1bSJed Brown /*TEST
154*c4762a1bSJed Brown 
155*c4762a1bSJed Brown    test:
156*c4762a1bSJed Brown 
157*c4762a1bSJed Brown TEST*/
158