xref: /petsc/src/sys/tests/ex40.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Test PETSc integer hash map.\n\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown #include <petsc/private/hashmapi.h>
4*c4762a1bSJed Brown #include <petsc/private/hashmapiv.h>
5*c4762a1bSJed Brown #include <petscsys.h>
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown /* Unused, keep it for testing purposes */
8*c4762a1bSJed Brown PETSC_HASH_MAP(HMapIP, PetscInt, void*, PetscHashInt, PetscHashEqual, NULL)
9*c4762a1bSJed Brown 
10*c4762a1bSJed Brown /* Unused, keep it for testing purposes */
11*c4762a1bSJed Brown typedef struct { double x; double y; double z; } Point;
12*c4762a1bSJed Brown static Point origin = {0.0, 0.0, 0.0};
13*c4762a1bSJed Brown PETSC_HASH_MAP(HMapIS, PetscInt, Point, PetscHashInt, PetscHashEqual, origin)
14*c4762a1bSJed Brown 
15*c4762a1bSJed Brown #define PetscAssert(expr) do {            \
16*c4762a1bSJed Brown if (PetscUnlikely(!(expr)))               \
17*c4762a1bSJed Brown   SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, \
18*c4762a1bSJed Brown            "Assertion: `%s' failed.",     \
19*c4762a1bSJed Brown            PetscStringize(expr));         \
20*c4762a1bSJed Brown } while(0)
21*c4762a1bSJed Brown 
22*c4762a1bSJed Brown 
23*c4762a1bSJed Brown int main(int argc,char **argv)
24*c4762a1bSJed Brown {
25*c4762a1bSJed Brown   PetscHMapI     ht = NULL, hd;
26*c4762a1bSJed Brown   PetscHMapIV    htv;
27*c4762a1bSJed Brown   PetscInt       n, v, koff, keys[4], voff, vals[4],na,nb,i,size,*karray,off;
28*c4762a1bSJed Brown   PetscScalar    *varray,*vwork;
29*c4762a1bSJed Brown   PetscBool      has, flag;
30*c4762a1bSJed Brown   PetscErrorCode ierr;
31*c4762a1bSJed Brown 
32*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
33*c4762a1bSJed Brown 
34*c4762a1bSJed Brown   ierr = PetscHMapICreate(&ht);CHKERRQ(ierr);
35*c4762a1bSJed Brown   PetscAssert(ht != NULL);
36*c4762a1bSJed Brown   ierr = PetscHMapIGetSize(ht,&n);CHKERRQ(ierr);
37*c4762a1bSJed Brown   PetscAssert(n == 0);
38*c4762a1bSJed Brown 
39*c4762a1bSJed Brown   ierr = PetscHMapIResize(ht,0);CHKERRQ(ierr);
40*c4762a1bSJed Brown   ierr = PetscHMapIGetSize(ht,&n);CHKERRQ(ierr);
41*c4762a1bSJed Brown   PetscAssert(n == 0);
42*c4762a1bSJed Brown 
43*c4762a1bSJed Brown   ierr = PetscHMapIHas(ht,123,&has);CHKERRQ(ierr);
44*c4762a1bSJed Brown   PetscAssert(has == PETSC_FALSE);
45*c4762a1bSJed Brown   ierr = PetscHMapIGet(ht,123,&v);CHKERRQ(ierr);
46*c4762a1bSJed Brown   PetscAssert(v == -1);
47*c4762a1bSJed Brown 
48*c4762a1bSJed Brown   ierr = PetscHMapISet(ht,123,42);CHKERRQ(ierr);
49*c4762a1bSJed Brown   ierr = PetscHMapIGetSize(ht,&n);CHKERRQ(ierr);
50*c4762a1bSJed Brown   PetscAssert(n == 1);
51*c4762a1bSJed Brown   ierr = PetscHMapIHas(ht,123,&has);CHKERRQ(ierr);
52*c4762a1bSJed Brown   PetscAssert(has == PETSC_TRUE);
53*c4762a1bSJed Brown   ierr = PetscHMapIGet(ht,123,&v);CHKERRQ(ierr);
54*c4762a1bSJed Brown   PetscAssert(v == 42);
55*c4762a1bSJed Brown 
56*c4762a1bSJed Brown   ierr = PetscHMapIDel(ht,123);CHKERRQ(ierr);
57*c4762a1bSJed Brown   ierr = PetscHMapIGetSize(ht,&n);CHKERRQ(ierr);
58*c4762a1bSJed Brown   PetscAssert(n == 0);
59*c4762a1bSJed Brown   ierr = PetscHMapIHas(ht,123,&has);CHKERRQ(ierr);
60*c4762a1bSJed Brown   PetscAssert(has == PETSC_FALSE);
61*c4762a1bSJed Brown   ierr = PetscHMapIGet(ht,123,&v);CHKERRQ(ierr);
62*c4762a1bSJed Brown   PetscAssert(v == -1);
63*c4762a1bSJed Brown 
64*c4762a1bSJed Brown   ierr = PetscHMapIQuerySet(ht,123,1,&flag);CHKERRQ(ierr);
65*c4762a1bSJed Brown   PetscAssert(flag == PETSC_TRUE);
66*c4762a1bSJed Brown   ierr = PetscHMapIQuerySet(ht,123,1,&flag);CHKERRQ(ierr);
67*c4762a1bSJed Brown   PetscAssert(flag == PETSC_FALSE);
68*c4762a1bSJed Brown   ierr = PetscHMapIQueryDel(ht,123,&flag);CHKERRQ(ierr);
69*c4762a1bSJed Brown   PetscAssert(flag == PETSC_TRUE);
70*c4762a1bSJed Brown   ierr = PetscHMapIQueryDel(ht,123,&flag);CHKERRQ(ierr);
71*c4762a1bSJed Brown   PetscAssert(flag == PETSC_FALSE);
72*c4762a1bSJed Brown 
73*c4762a1bSJed Brown   ierr = PetscHMapIResize(ht,13);CHKERRQ(ierr);
74*c4762a1bSJed Brown   ierr = PetscHMapIGetSize(ht,&n);CHKERRQ(ierr);
75*c4762a1bSJed Brown   PetscAssert(n == 0);
76*c4762a1bSJed Brown 
77*c4762a1bSJed Brown   ierr = PetscHMapIClear(ht);CHKERRQ(ierr);
78*c4762a1bSJed Brown   ierr = PetscHMapIGetSize(ht,&n);CHKERRQ(ierr);
79*c4762a1bSJed Brown   PetscAssert(n == 0);
80*c4762a1bSJed Brown 
81*c4762a1bSJed Brown   ierr = PetscHMapISet(ht,321,24);CHKERRQ(ierr);
82*c4762a1bSJed Brown   ierr = PetscHMapISet(ht,123,42);CHKERRQ(ierr);
83*c4762a1bSJed Brown   ierr = PetscHMapIGetSize(ht,&n);CHKERRQ(ierr);
84*c4762a1bSJed Brown   PetscAssert(n == 2);
85*c4762a1bSJed Brown 
86*c4762a1bSJed Brown   koff = 0; keys[0] = keys[1] = 0;CHKERRQ(ierr);
87*c4762a1bSJed Brown   ierr = PetscHMapIGetKeys(ht,&koff,keys);CHKERRQ(ierr);
88*c4762a1bSJed Brown   ierr = PetscSortInt(koff,keys);CHKERRQ(ierr);
89*c4762a1bSJed Brown   PetscAssert(koff == 2);
90*c4762a1bSJed Brown   PetscAssert(keys[0] == 123);
91*c4762a1bSJed Brown   PetscAssert(keys[1] == 321);
92*c4762a1bSJed Brown 
93*c4762a1bSJed Brown   voff = 0; vals[0] = vals[1] = 0;
94*c4762a1bSJed Brown   ierr = PetscHMapIGetVals(ht,&voff,vals);CHKERRQ(ierr);
95*c4762a1bSJed Brown   ierr = PetscSortInt(voff,vals);CHKERRQ(ierr);
96*c4762a1bSJed Brown   PetscAssert(voff == 2);
97*c4762a1bSJed Brown   PetscAssert(vals[0] == 24);
98*c4762a1bSJed Brown   PetscAssert(vals[1] == 42);
99*c4762a1bSJed Brown 
100*c4762a1bSJed Brown   koff = 0; keys[0] = keys[1] = 0;
101*c4762a1bSJed Brown   voff = 0; vals[0] = vals[1] = 0;
102*c4762a1bSJed Brown   ierr = PetscHMapIDuplicate(ht,&hd);CHKERRQ(ierr);
103*c4762a1bSJed Brown   ierr = PetscHMapIGetKeys(ht,&koff,keys);CHKERRQ(ierr);
104*c4762a1bSJed Brown   ierr = PetscHMapIGetVals(ht,&voff,vals);CHKERRQ(ierr);
105*c4762a1bSJed Brown   ierr = PetscSortInt(koff,keys);CHKERRQ(ierr);
106*c4762a1bSJed Brown   ierr = PetscSortInt(voff,vals);CHKERRQ(ierr);
107*c4762a1bSJed Brown   PetscAssert(koff == 2);
108*c4762a1bSJed Brown   PetscAssert(voff == 2);
109*c4762a1bSJed Brown   PetscAssert(keys[0] == 123);
110*c4762a1bSJed Brown   PetscAssert(keys[1] == 321);
111*c4762a1bSJed Brown   PetscAssert(vals[0] == 24);
112*c4762a1bSJed Brown   PetscAssert(vals[1] == 42);
113*c4762a1bSJed Brown   ierr = PetscHMapIDestroy(&hd);CHKERRQ(ierr);
114*c4762a1bSJed Brown 
115*c4762a1bSJed Brown   ierr = PetscHMapISet(ht,0,0);CHKERRQ(ierr);
116*c4762a1bSJed Brown   ierr = PetscHMapIGetSize(ht,&n);CHKERRQ(ierr);
117*c4762a1bSJed Brown   PetscAssert(n != 0);
118*c4762a1bSJed Brown   ierr = PetscHMapIReset(ht);CHKERRQ(ierr);
119*c4762a1bSJed Brown   ierr = PetscHMapIGetSize(ht,&n);CHKERRQ(ierr);
120*c4762a1bSJed Brown   PetscAssert(n == 0);
121*c4762a1bSJed Brown   ierr = PetscHMapIReset(ht);CHKERRQ(ierr);
122*c4762a1bSJed Brown   ierr = PetscHMapIGetSize(ht,&n);CHKERRQ(ierr);
123*c4762a1bSJed Brown   PetscAssert(n == 0);
124*c4762a1bSJed Brown   ierr = PetscHMapISet(ht,0,0);CHKERRQ(ierr);
125*c4762a1bSJed Brown   ierr = PetscHMapIGetSize(ht,&n);CHKERRQ(ierr);
126*c4762a1bSJed Brown   PetscAssert(n != 0);
127*c4762a1bSJed Brown 
128*c4762a1bSJed Brown   ierr = PetscHMapIDestroy(&ht);CHKERRQ(ierr);
129*c4762a1bSJed Brown   PetscAssert(ht == NULL);
130*c4762a1bSJed Brown 
131*c4762a1bSJed Brown   ierr = PetscHMapICreate(&ht);CHKERRQ(ierr);
132*c4762a1bSJed Brown   ierr = PetscHMapIReset(ht);CHKERRQ(ierr);
133*c4762a1bSJed Brown   ierr = PetscHMapIGetSize(ht,&n);CHKERRQ(ierr);
134*c4762a1bSJed Brown   PetscAssert(n == 0);
135*c4762a1bSJed Brown   ierr = PetscHMapIDestroy(&ht);CHKERRQ(ierr);
136*c4762a1bSJed Brown 
137*c4762a1bSJed Brown   ierr = PetscHMapIVCreate(&htv);CHKERRQ(ierr);
138*c4762a1bSJed Brown   n = 10;
139*c4762a1bSJed Brown   ierr = PetscHMapIVResize(htv,n);CHKERRQ(ierr);
140*c4762a1bSJed Brown   ierr = PetscHMapIVGetCapacity(htv,&na);CHKERRQ(ierr);
141*c4762a1bSJed Brown   PetscAssert(na>=n);
142*c4762a1bSJed Brown   for (i=0; i<n; i++) {
143*c4762a1bSJed Brown     ierr = PetscHMapIVSet(htv,i+100,10.);CHKERRQ(ierr);
144*c4762a1bSJed Brown   }
145*c4762a1bSJed Brown   ierr = PetscHMapIVGetCapacity(htv,&nb);CHKERRQ(ierr);
146*c4762a1bSJed Brown   PetscAssert(nb>=na);
147*c4762a1bSJed Brown   for (i=0; i<(2*n); i++) {
148*c4762a1bSJed Brown     ierr = PetscHMapIVAddValue(htv,i+100,5.);CHKERRQ(ierr);
149*c4762a1bSJed Brown   }
150*c4762a1bSJed Brown   ierr = PetscHMapIVGetSize(htv,&size);CHKERRQ(ierr);
151*c4762a1bSJed Brown   PetscAssert(size==(2*n));
152*c4762a1bSJed Brown   ierr = PetscMalloc3(size,&karray,size,&varray,size,&vwork);CHKERRQ(ierr);
153*c4762a1bSJed Brown   off = 0;
154*c4762a1bSJed Brown   ierr = PetscHMapIVGetPairs(htv,&off,karray,varray);CHKERRQ(ierr);
155*c4762a1bSJed Brown   PetscAssert(off==(2*n));
156*c4762a1bSJed Brown   ierr = PetscSortIntWithDataArray(off,karray,varray,sizeof(PetscScalar),vwork);CHKERRQ(ierr);
157*c4762a1bSJed Brown   for (i=0; i<n; i++) {
158*c4762a1bSJed Brown     PetscAssert(karray[i]==(i+100));
159*c4762a1bSJed Brown     PetscAssert(karray[n+i]==(n+i+100));
160*c4762a1bSJed Brown     PetscAssert(varray[i]==15.);
161*c4762a1bSJed Brown     PetscAssert(varray[n+i]==5.);
162*c4762a1bSJed Brown   }
163*c4762a1bSJed Brown   ierr = PetscFree3(karray,varray,vwork);CHKERRQ(ierr);
164*c4762a1bSJed Brown   ierr = PetscHMapIVDestroy(&htv);CHKERRQ(ierr);
165*c4762a1bSJed Brown 
166*c4762a1bSJed Brown   ierr = PetscFinalize();
167*c4762a1bSJed Brown   return ierr;
168*c4762a1bSJed Brown }
169*c4762a1bSJed Brown 
170*c4762a1bSJed Brown 
171*c4762a1bSJed Brown /*TEST
172*c4762a1bSJed Brown 
173*c4762a1bSJed Brown    test:
174*c4762a1bSJed Brown 
175*c4762a1bSJed Brown TEST*/
176