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