Actual source code: hashtable.h
1: #if !defined(PETSC_HASHTABLE_H)
2: #define PETSC_HASHTABLE_H
4: #include <petsc/private/petscimpl.h>
6: #define kh_inline inline
7: #define klib_unused PETSC_UNUSED
8: #include <petsc/private/khash/khash.h>
10: /* Required for khash <= 0.2.5 */
11: #if !defined(kcalloc)
12: #define kcalloc(N,Z) calloc(N,Z)
13: #endif
14: #if !defined(kmalloc)
15: #define kmalloc(Z) malloc(Z)
16: #endif
17: #if !defined(krealloc)
18: #define krealloc(P,Z) realloc(P,Z)
19: #endif
20: #if !defined(kfree)
21: #define kfree(P) free(P)
22: #endif
24: /* --- Useful extensions to khash --- */
26: #if !defined(kh_reset)
27: /*! @function
28: @abstract Reset a hash table to initial state.
29: @param name Name of the hash table [symbol]
30: @param h Pointer to the hash table [khash_t(name)*]
31: */
32: #define kh_reset(name, h) { \
33: if (h) { \
34: kfree((h)->keys); kfree((h)->flags); \
35: kfree((h)->vals); \
36: memset((h), 0x00, sizeof(*(h))); \
37: } }
38: #endif /*kh_reset*/
40: #if !defined(kh_foreach)
41: /*! @function
42: @abstract Iterate over the entries in the hash table
43: @param h Pointer to the hash table [khash_t(name)*]
44: @param kvar Variable to which key will be assigned
45: @param vvar Variable to which value will be assigned
46: @param code Block of code to execute
47: */
48: #define kh_foreach(h, kvar, vvar, code) { khint_t __i; \
49: for (__i = kh_begin(h); __i != kh_end(h); ++__i) { \
50: if (!kh_exist(h,__i)) continue; \
51: (kvar) = kh_key(h,__i); \
52: (vvar) = kh_val(h,__i); \
53: code; \
54: } }
55: #endif /*kh_foreach*/
57: #if !defined(kh_foreach_key)
58: /*! @function
59: @abstract Iterate over the keys in the hash table
60: @param h Pointer to the hash table [khash_t(name)*]
61: @param kvar Variable to which key will be assigned
62: @param code Block of code to execute
63: */
64: #define kh_foreach_key(h, kvar, code) { khint_t __i; \
65: for (__i = kh_begin(h); __i != kh_end(h); ++__i) { \
66: if (!kh_exist(h,__i)) continue; \
67: (kvar) = kh_key(h,__i); \
68: code; \
69: } }
70: #endif /*kh_foreach_key*/
72: #if !defined(kh_foreach_value)
73: /*! @function
74: @abstract Iterate over the values in the hash table
75: @param h Pointer to the hash table [khash_t(name)*]
76: @param vvar Variable to which value will be assigned
77: @param code Block of code to execute
78: */
79: #define kh_foreach_value(h, vvar, code) { khint_t __i; \
80: for (__i = kh_begin(h); __i != kh_end(h); ++__i) { \
81: if (!kh_exist(h,__i)) continue; \
82: (vvar) = kh_val(h,__i); \
83: code; \
84: } }
85: #endif /*kh_foreach_value*/
87: /* --- Helper macro for error checking --- */
89: #if defined(PETSC_USE_DEBUG)
91: #else
92: #define PetscHashAssert(expr) ((void)(expr))
93: #endif
95: /* --- Low level iterator API --- */
97: typedef khiter_t PetscHashIter;
99: #define PetscHashIterBegin(ht,i) do { \
100: (i) = kh_begin((ht)); \
101: if ((i) != kh_end((ht)) && !kh_exist((ht),(i))) \
102: PetscHashIterNext((ht),(i)); \
103: } while (0)
105: #define PetscHashIterNext(ht,i) \
106: do { ++(i); } while ((i) != kh_end((ht)) && !kh_exist((ht),(i)))
108: #define PetscHashIterAtEnd(ht,i) ((i) == kh_end((ht)))
110: #define PetscHashIterGetKey(ht,i,k) ((k) = kh_key((ht),(i)))
112: #define PetscHashIterGetVal(ht,i,v) ((v) = kh_val((ht),(i)))
114: #define PetscHashIterSetVal(ht,i,v) (kh_val((ht),(i)) = (v))
116: /* --- Thomas Wang integer hash functions --- */
118: typedef khint32_t PetscHash32_t;
119: typedef khint64_t PetscHash64_t;
120: typedef khint_t PetscHash_t;
122: /* Thomas Wang's first version for 32bit integers */
123: static inline PetscHash_t PetscHash_UInt32_v0(PetscHash32_t key)
124: {
125: key += ~(key << 15);
126: key ^= (key >> 10);
127: key += (key << 3);
128: key ^= (key >> 6);
129: key += ~(key << 11);
130: key ^= (key >> 16);
131: return key;
132: }
134: /* Thomas Wang's second version for 32bit integers */
135: static inline PetscHash_t PetscHash_UInt32_v1(PetscHash32_t key)
136: {
137: key = ~key + (key << 15); /* key = (key << 15) - key - 1; */
138: key = key ^ (key >> 12);
139: key = key + (key << 2);
140: key = key ^ (key >> 4);
141: key = key * 2057; /* key = (key + (key << 3)) + (key << 11); */
142: key = key ^ (key >> 16);
143: return key;
144: }
146: static inline PetscHash_t PetscHash_UInt32(PetscHash32_t key)
147: {
148: return PetscHash_UInt32_v1(key);
149: }
151: /* Thomas Wang's version for 64bit integer -> 32bit hash */
152: static inline PetscHash32_t PetscHash_UInt64_32(PetscHash64_t key)
153: {
154: key = ~key + (key << 18); /* key = (key << 18) - key - 1; */
155: key = key ^ (key >> 31);
156: key = key * 21; /* key = (key + (key << 2)) + (key << 4); */
157: key = key ^ (key >> 11);
158: key = key + (key << 6);
159: key = key ^ (key >> 22);
160: return (PetscHash32_t)key;
161: }
163: /* Thomas Wang's version for 64bit integer -> 64bit hash */
164: static inline PetscHash64_t PetscHash_UInt64_64(PetscHash64_t key)
165: {
166: key = ~key + (key << 21); /* key = (key << 21) - key - 1; */
167: key = key ^ (key >> 24);
168: key = key * 265; /* key = (key + (key << 3)) + (key << 8); */
169: key = key ^ (key >> 14);
170: key = key * 21; /* key = (key + (key << 2)) + (key << 4); */
171: key = key ^ (key >> 28);
172: key = key + (key << 31);
173: return key;
174: }
176: static inline PetscHash_t PetscHash_UInt64(PetscHash64_t key)
177: {
178: return sizeof(PetscHash_t) < sizeof(PetscHash64_t)
179: ? (PetscHash_t)PetscHash_UInt64_32(key)
180: : (PetscHash_t)PetscHash_UInt64_64(key);
181: }
183: static inline PetscHash_t PetscHashInt(PetscInt key)
184: {
185: #if defined(PETSC_USE_64BIT_INDICES)
186: return PetscHash_UInt64((PetscHash64_t)key);
187: #else
188: return PetscHash_UInt32((PetscHash32_t)key);
189: #endif
190: }
192: static inline PetscHash_t PetscHashPointer(void *key)
193: {
194: #if PETSC_SIZEOF_VOID_P == 8
195: return PetscHash_UInt64((PetscHash64_t)key);
196: #else
197: return PetscHash_UInt32((PetscHash32_t)key);
198: #endif
199: }
201: static inline PetscHash_t PetscHashCombine(PetscHash_t seed, PetscHash_t hash)
202: {
203: /* https://doi.org/10.1002/asi.10170 */
204: /* https://dl.acm.org/citation.cfm?id=759509 */
205: return seed ^ (hash + (seed << 6) + (seed >> 2));
206: }
208: #define PetscHashEqual(a,b) ((a) == (b))
210: #endif /* PETSC_HASHTABLE_H */