Hasher.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Description
27  Hashing functions, mostly from Bob Jenkins
28 \*---------------------------------------------------------------------------*/
29 
30 #include "Hasher.H"
31 #include "HasherInt.H"
32 #include "endian.H"
33 
34 // Left-rotate a 32-bit value and carry by nBits
35 #define bitRotateLeft(x, nBits) (((x) << (nBits)) | ((x) >> (32 - (nBits))))
36 
37 
38 // ----------------------------------------------------------------------------
39 // lookup3.c, by Bob Jenkins, May 2006, Public Domain.
40 //
41 // These are functions for producing 32-bit hashes for hash table lookup.
42 // hashword(), hashlittle(), hashlittle2(), hashbig(), mix(), and final()
43 // are externally useful functions. Routines to test the hash are included
44 // if SELF_TEST is defined. You can use this free for any purpose. It's in
45 // the public domain. It has no warranty.
46 //
47 // You probably want to use hashlittle(). hashlittle() and hashbig()
48 // hash byte arrays. hashlittle() is is faster than hashbig() on
49 // little-endian machines. Intel and AMD are little-endian machines.
50 // On second thought, you probably want hashlittle2(), which is identical to
51 // hashlittle() except it returns two 32-bit hashes for the price of one.
52 // You could implement hashbig2() if you wanted but I haven't bothered here.
53 //
54 // If you want to find a hash of, say, exactly 7 integers, do
55 // a = i1; b = i2; c = i3;
56 // mix(a,b,c);
57 // a += i4; b += i5; c += i6;
58 // mix(a,b,c);
59 // a += i7;
60 // final(a,b,c);
61 // then use c as the hash value. If you have a variable length array of
62 // 4-byte integers to hash, use hashword(). If you have a byte array (like
63 // a character string), use hashlittle(). If you have several byte arrays, or
64 // a mix of things, see the comments above hashlittle().
65 //
66 // Why is this so big? I read 12 bytes at a time into 3 4-byte integers,
67 // then mix those integers. This is fast (you can do a lot more thorough
68 // mixing with 12*3 instructions on 3 integers than you can with 3 instructions
69 // on 1 byte), but shoehorning those bytes into integers efficiently is messy.
70 // ----------------------------------------------------------------------------
71 
72 // ----------------------------------------------------------------------------
73 // mix -- mix 3 32-bit values reversibly.
74 //
75 // This is reversible, so any information in (a,b,c) before mix_hash() is
76 // still in (a,b,c) after mix_hash().
77 //
78 // If four pairs of (a,b,c) inputs are run through mix_hash(), or through
79 // mix_hash() in reverse, there are at least 32 bits of the output that
80 // are sometimes the same for one pair and different for another pair.
81 // This was tested for:
82 // * pairs that differed by one bit, by two bits, in any combination
83 // of top bits of (a,b,c), or in any combination of bottom bits of
84 // (a,b,c).
85 // * "differ" is defined as +, -, ^, or ~^. For + and -, I transformed
86 // the output delta to a Gray code (a^(a>>1)) so a string of 1's (as
87 // is commonly produced by subtraction) look like a single 1-bit
88 // difference.
89 // * the base values were pseudorandom, all zero but one bit set, or
90 // all zero plus a counter that starts at zero.
91 //
92 // Some k values for my "a-=c; a^=rot(c,k); c+=b;" arrangement that
93 // satisfy this are
94 // 4 6 8 16 19 4
95 // 9 15 3 18 27 15
96 // 14 9 3 7 17 3
97 // Well, "9 15 3 18 27 15" didn't quite get 32 bits diffing
98 // for "differ" defined as + with a one-bit base and a two-bit delta. I
99 // used http://burtleburtle.net/bob/hash/avalanche.html to choose
100 // the operations, constants, and arrangements of the variables.
101 //
102 // This does not achieve avalanche. There are input bits of (a,b,c)
103 // that fail to affect some output bits of (a,b,c), especially of a. The
104 // most thoroughly mixed value is c, but it doesn't really even achieve
105 // avalanche in c.
106 //
107 // This allows some parallelism. Read-after-writes are good at doubling
108 // the number of bits affected, so the goal of mixing pulls in the opposite
109 // direction as the goal of parallelism. I did what I could. Rotates
110 // seem to cost as much as shifts on every machine I could lay my hands
111 // on, and rotates are much kinder to the top and bottom bits, so I used
112 // rotates.
113 // ----------------------------------------------------------------------------
114 
115 #define bitMixer(a, b, c) \
116  { \
117  a -= c; a ^= bitRotateLeft(c, 4); c += b; \
118  b -= a; b ^= bitRotateLeft(a, 6); a += c; \
119  c -= b; c ^= bitRotateLeft(b, 8); b += a; \
120  a -= c; a ^= bitRotateLeft(c,16); c += b; \
121  b -= a; b ^= bitRotateLeft(a,19); a += c; \
122  c -= b; c ^= bitRotateLeft(b, 4); b += a; \
123  }
124 
125 
126 // ----------------------------------------------------------------------------
127 // final -- final mixing of 3 32-bit values (a,b,c) into c
128 //
129 // Pairs of (a,b,c) values differing in only a few bits will usually
130 // produce values of c that look totally different. This was tested for
131 // * pairs that differed by one bit, by two bits, in any combination
132 // of top bits of (a,b,c), or in any combination of bottom bits of
133 // (a,b,c).
134 // * "differ" is defined as +, -, ^, or ~^. For + and -, I transformed
135 // the output delta to a Gray code (a^(a>>1)) so a string of 1's (as
136 // is commonly produced by subtraction) look like a single 1-bit
137 // difference.
138 // * the base values were pseudorandom, all zero but one bit set, or
139 // all zero plus a counter that starts at zero.
140 //
141 // These constants passed:
142 // 14 11 25 16 4 14 24
143 // 12 14 25 16 4 14 24
144 // and these came close:
145 // 4 8 15 26 3 22 24
146 // 10 8 15 26 3 22 24
147 // 11 8 15 26 3 22 24
148 // ----------------------------------------------------------------------------
149 
150 #define bitMixerFinal(a, b, c) \
151  { \
152  c ^= b; c -= bitRotateLeft(b, 14); \
153  a ^= c; a -= bitRotateLeft(c, 11); \
154  b ^= a; b -= bitRotateLeft(a, 25); \
155  c ^= b; c -= bitRotateLeft(b, 16); \
156  a ^= c; a -= bitRotateLeft(c, 4); \
157  b ^= a; b -= bitRotateLeft(a, 14); \
158  c ^= b; c -= bitRotateLeft(b, 24); \
159  }
160 
161 
162 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
163 
164 // ----------------------------------------------------------------------------
165 // This works on all machines. To be useful, it requires
166 // -- that the key be an array of uint32_t's, and
167 // -- that the length be the number of uint32_t's in the key
168 //
169 // The function hashword() is identical to hashlittle() on little-endian
170 // machines, and identical to hashbig() on big-endian machines,
171 // except that the length has to be measured in uint32_ts rather than in
172 // bytes. hashlittle() is more complicated than hashword() only because
173 // hashlittle() has to dance around fitting the key bytes into registers.
174 // ----------------------------------------------------------------------------
175 unsigned Foam::HasherInt
176 (
177  const uint32_t *k,
178  size_t length,
179  unsigned seed
180 )
181 {
182  uint32_t a, b, c;
183 
184  // Set up the internal state
185  a = b = c = 0xdeadbeef + (static_cast<uint32_t>(length) << 2) + seed;
186 
187  // handle most of the key
188  while (length > 3)
189  {
190  a += k[0];
191  b += k[1];
192  c += k[2];
193  bitMixer(a,b,c);
194  length -= 3;
195  k += 3;
196  }
197 
198  // handle the last 3 uint32_t's
199  switch (length) // all case statements fall through
200  {
201  case 3 : c += k[2]; [[fallthrough]];
202  case 2 : b += k[1]; [[fallthrough]];
203  case 1 : a += k[0];
204  bitMixerFinal(a,b,c);
205  [[fallthrough]];
206  case 0 : // case 0: nothing left to add
207  break;
208  }
209 
210  return c;
211 }
212 
213 
214 // ----------------------------------------------------------------------------
215 // hashword2() -- same as hashword(), but take two seeds and return two
216 // 32-bit values. pc and pb must both be non-null, and *pc and *pb must
217 // both be initialized with seeds. If you pass in (*pb)==0, the output
218 // (*pc) will be the same as the return value from hashword().
219 // ----------------------------------------------------------------------------
220 
221 // Currently unused
222 #if 0
223 unsigned Foam::HasherDual
224 (
225  const uint32_t *k,
226  size_t length,
227  unsigned& hash1, // IN: seed OUT: primary hash value
228  unsigned& hash2 // IN: more seed OUT: secondary hash value
229 )
230 {
231  uint32_t a, b, c;
232 
233  // Set up the internal state
234  a = b = c = 0xdeadbeef + (static_cast<uint32_t>(length) << 2) + hash1;
235  c += hash2;
236 
237  // handle most of the key
238  while (length > 3)
239  {
240  a += k[0];
241  b += k[1];
242  c += k[2];
243  bitMixer(a,b,c);
244  length -= 3;
245  k += 3;
246  }
247 
248  // handle the last 3 uint32_t's
249  switch (length) // all case statements fall through
250  {
251  case 3 : c += k[2]; [[fallthrough]];
252  case 2 : b += k[1]; [[fallthrough]];
253  case 1 : a += k[0];
254  bitMixerFinal(a,b,c);
255  [[fallthrough]];
256  case 0 : // case 0: nothing left to add
257  break;
258  }
259 
260  // report the result
261  hash1 = c;
262  hash2 = b;
263 
264  // return primary hash value
265  return c;
266 }
267 #endif
268 
269 
270 // * * * * * * * * * * * * * * Static Functions * * * * * * * * * * * * * * //
271 
272 // ----------------------------------------------------------------------------
273 // hashlittle() -- hash a variable-length key into a 32-bit value
274 // k : the key (the unaligned variable-length array of bytes)
275 // length : the length of the key, counting by bytes
276 // initval : can be any 4-byte value
277 // Returns a 32-bit value. Every bit of the key affects every bit of
278 // the return value. Two keys differing by one or two bits will have
279 // totally different hash values.
280 //
281 // The best hash table sizes are powers of 2. There is no need to do
282 // mod a prime (mod is sooo slow!). If you need less than 32 bits,
283 // use a bitmask. For example, if you need only 10 bits, do
284 // h = (h & hashmask(10));
285 // In which case, the hash table should have hashsize(10) elements.
286 //
287 // If you are hashing n strings (uint8_t **)k, do it like this:
288 // for (i=0, h=0; i<n; ++i) h = hashlittle( k[i], len[i], h);
289 //
290 // By Bob Jenkins, 2006. bob_jenkins@burtleburtle.net. You may use this
291 // code any way you wish, private, educational, or commercial. It's free.
292 //
293 // Use for hash table lookup, or anything where one collision in 2^^32 is
294 // acceptable. Do NOT use for cryptographic purposes.
295 // ----------------------------------------------------------------------------
296 
297 // Specialized little-endian code
298 #ifdef WM_LITTLE_ENDIAN
299 static unsigned jenkins_hashlittle
300 (
301  const void *key,
302  size_t length,
303  unsigned initval
304 )
305 {
306  uint32_t a, b, c;
307  union { const void *ptr; size_t i; } u; // to cast key to (size_t) happily
308 
309  // Set up the internal state
310  a = b = c = 0xdeadbeef + static_cast<uint32_t>(length) + initval;
311 
312  u.ptr = key;
313  if ((u.i & 0x3) == 0)
314  {
315  // 32-bit chunks
316  const uint32_t *k = reinterpret_cast<const uint32_t*>(key);
317 
318  // all but last block: aligned reads and affect 32 bits of (a,b,c)
319  while (length > 12)
320  {
321  a += k[0];
322  b += k[1];
323  c += k[2];
324  bitMixer(a,b,c);
325  length -= 12;
326  k += 3;
327  }
328 
329  // handle the last (probably partial) block byte-wise
330  const uint8_t *k8 = reinterpret_cast<const uint8_t*>(k);
331  switch (length)
332  {
333  case 12: c += k[2]; b += k[1]; a += k[0]; break;
334  case 11: c += static_cast<uint32_t>(k8[10]) << 16; [[fallthrough]];
335  case 10: c += static_cast<uint32_t>(k8[9]) << 8; [[fallthrough]];
336  case 9 : c += k8[8]; [[fallthrough]];
337  case 8 : b += k[1]; a += k[0]; break;
338  case 7 : b += static_cast<uint32_t>(k8[6]) << 16; [[fallthrough]];
339  case 6 : b += static_cast<uint32_t>(k8[5]) << 8; [[fallthrough]];
340  case 5 : b += k8[4]; [[fallthrough]];
341  case 4 : a += k[0]; break;
342  case 3 : a += static_cast<uint32_t>(k8[2]) << 16; [[fallthrough]];
343  case 2 : a += static_cast<uint32_t>(k8[1]) << 8; [[fallthrough]];
344  case 1 : a += k8[0]; break;
345  case 0 : return c; // zero-length requires no mixing
346  }
347  }
348  else if ((u.i & 0x1) == 0)
349  {
350  // 16-bit chunks
351  const uint16_t *k = reinterpret_cast<const uint16_t*>(key);
352 
353  // all but last block: aligned reads and different mixing
354  while (length > 12)
355  {
356  a += k[0] + (static_cast<uint32_t>(k[1]) << 16);
357  b += k[2] + (static_cast<uint32_t>(k[3]) << 16);
358  c += k[4] + (static_cast<uint32_t>(k[5]) << 16);
359  bitMixer(a,b,c);
360  length -= 12;
361  k += 6;
362  }
363 
364  // handle the last (probably partial) block
365  const uint8_t *k8 = reinterpret_cast<const uint8_t*>(k);
366  switch (length)
367  {
368  case 12:
369  c += k[4] + (static_cast<uint32_t>(k[5]) << 16);
370  b += k[2] + (static_cast<uint32_t>(k[3]) << 16);
371  a += k[0] + (static_cast<uint32_t>(k[1]) << 16);
372  break;
373  case 11:
374  c += static_cast<uint32_t>(k8[10]) << 16;
375  [[fallthrough]];
376  case 10:
377  c += k[4];
378  b += k[2] + (static_cast<uint32_t>(k[3]) << 16);
379  a += k[0] + (static_cast<uint32_t>(k[1]) << 16);
380  break;
381  case 9 :
382  c += k8[8];
383  [[fallthrough]];
384  case 8 :
385  b += k[2] + (static_cast<uint32_t>(k[3]) << 16);
386  a += k[0] + (static_cast<uint32_t>(k[1]) << 16);
387  break;
388  case 7 :
389  b += static_cast<uint32_t>(k8[6]) << 16;
390  [[fallthrough]];
391  case 6 :
392  b += k[2];
393  a += k[0] + (static_cast<uint32_t>(k[1]) << 16);
394  break;
395  case 5 :
396  b += k8[4];
397  [[fallthrough]];
398  case 4 :
399  a += k[0] + (static_cast<uint32_t>(k[1]) << 16);
400  break;
401  case 3 :
402  a += static_cast<uint32_t>(k8[2]) << 16;
403  [[fallthrough]];
404  case 2 :
405  a += k[0];
406  break;
407  case 1 :
408  a += k8[0];
409  break;
410  case 0 : return c; // zero-length requires no mixing
411  }
412  }
413  else
414  {
415  // need to read the key one byte at a time
416  const uint8_t *k = reinterpret_cast<const uint8_t*>(key);
417 
418  // all but the last block: affect some 32 bits of (a,b,c)
419  while (length > 12)
420  {
421  a += k[0];
422  a += static_cast<uint32_t>(k[1]) << 8;
423  a += static_cast<uint32_t>(k[2]) << 16;
424  a += static_cast<uint32_t>(k[3]) << 24;
425  b += k[4];
426  b += static_cast<uint32_t>(k[5]) << 8;
427  b += static_cast<uint32_t>(k[6]) << 16;
428  b += static_cast<uint32_t>(k[7]) << 24;
429  c += k[8];
430  c += static_cast<uint32_t>(k[9]) << 8;
431  c += static_cast<uint32_t>(k[10]) << 16;
432  c += static_cast<uint32_t>(k[11]) << 24;
433 
434  bitMixer(a,b,c);
435  length -= 12;
436  k += 12;
437  }
438 
439  // last block: affect all 32 bits of (c)
440  switch (length) // most case statements fall through
441  {
442  case 12: c += static_cast<uint32_t>(k[11]) << 24; [[fallthrough]];
443  case 11: c += static_cast<uint32_t>(k[10]) << 16; [[fallthrough]];
444  case 10: c += static_cast<uint32_t>(k[9]) << 8; [[fallthrough]];
445  case 9 : c += k[8]; [[fallthrough]];
446 
447  case 8 : b += static_cast<uint32_t>(k[7]) << 24; [[fallthrough]];
448  case 7 : b += static_cast<uint32_t>(k[6]) << 16; [[fallthrough]];
449  case 6 : b += static_cast<uint32_t>(k[5]) << 8; [[fallthrough]];
450  case 5 : b += k[4]; [[fallthrough]];
451 
452  case 4 : a += static_cast<uint32_t>(k[3]) << 24; [[fallthrough]];
453  case 3 : a += static_cast<uint32_t>(k[2]) << 16; [[fallthrough]];
454  case 2 : a += static_cast<uint32_t>(k[1]) << 8; [[fallthrough]];
455  case 1 : a += k[0];
456  break;
457 
458  case 0 : return c;
459  }
460  }
461 
462  bitMixerFinal(a,b,c);
463  return c;
464 }
465 #endif
466 
467 
468 // ----------------------------------------------------------------------------
469 // hashbig():
470 // This is the same as hashword() on big-endian machines. It is different
471 // from hashlittle() on all machines. hashbig() takes advantage of
472 // big-endian byte ordering.
473 // ----------------------------------------------------------------------------
474 // specialized big-endian code
475 #ifdef WM_BIG_ENDIAN
476 static unsigned jenkins_hashbig
477 (
478  const void *key,
479  size_t length,
480  unsigned initval
481 )
482 {
483  uint32_t a, b, c;
484  union { const void *ptr; size_t i; } u; // to cast key to (size_t) happily
485 
486  // Set up the internal state
487  a = b = c = 0xdeadbeef + static_cast<uint32_t>(length) + initval;
488 
489  u.ptr = key;
490  if ((u.i & 0x3) == 0)
491  {
492  // 32-bit chunks
493  const uint32_t *k = reinterpret_cast<const uint32_t*>(key);
494 
495  // all but last block: aligned reads and affect 32 bits of (a,b,c)
496  while (length > 12)
497  {
498  a += k[0];
499  b += k[1];
500  c += k[2];
501  bitMixer(a,b,c);
502  length -= 12;
503  k += 3;
504  }
505 
506  // handle the last (probably partial) block byte-wise
507  const uint8_t *k8 = reinterpret_cast<const uint8_t*>(k);
508 
509  switch (length) // most the case statements fall through
510  {
511  case 12: c += k[2]; b += k[1]; a += k[0]; break;
512  case 11: c += static_cast<uint32_t>(k8[10]) << 8; [[fallthrough]];
513  case 10: c += static_cast<uint32_t>(k8[9]) << 16; [[fallthrough]];
514  case 9 : c += static_cast<uint32_t>(k8[8]) << 24; [[fallthrough]];
515  case 8 : b += k[1]; a += k[0]; break;
516  case 7 : b += static_cast<uint32_t>(k8[6]) << 8; [[fallthrough]];
517  case 6 : b += static_cast<uint32_t>(k8[5]) << 16; [[fallthrough]];
518  case 5 : b += static_cast<uint32_t>(k8[4]) << 24; [[fallthrough]];
519  case 4 : a += k[0]; break;
520  case 3 : a += static_cast<uint32_t>(k8[2]) << 8; [[fallthrough]];
521  case 2 : a += static_cast<uint32_t>(k8[1]) << 16; [[fallthrough]];
522  case 1 : a += static_cast<uint32_t>(k8[0]) << 24; break;
523  case 0 : return c;
524  }
525  }
526  else
527  {
528  // need to read the key one byte at a time
529  const uint8_t *k = reinterpret_cast<const uint8_t*>(key);
530 
531  // all but the last block: affect some 32 bits of (a,b,c)
532  while (length > 12)
533  {
534  a += static_cast<uint32_t>(k[0]) << 24;
535  a += static_cast<uint32_t>(k[1]) << 16;
536  a += static_cast<uint32_t>(k[2]) << 8;
537  a += static_cast<uint32_t>(k[3]);
538  b += static_cast<uint32_t>(k[4]) << 24;
539  b += static_cast<uint32_t>(k[5]) << 16;
540  b += static_cast<uint32_t>(k[6]) << 8;
541  b += static_cast<uint32_t>(k[7]);
542  c += static_cast<uint32_t>(k[8]) << 24;
543  c += static_cast<uint32_t>(k[9]) << 16;
544  c += static_cast<uint32_t>(k[10]) << 8;
545  c += static_cast<uint32_t>(k[11]);
546 
547  bitMixer(a,b,c);
548  length -= 12;
549  k += 12;
550  }
551 
552  // last block: affect all 32 bits of (c)
553  switch (length) // the case statements fall through
554  {
555  case 12: c += k[11]; [[fallthrough]];
556  case 11: c += static_cast<uint32_t>(k[10]) << 8; [[fallthrough]];
557  case 10: c += static_cast<uint32_t>(k[9]) << 16; [[fallthrough]];
558  case 9 : c += static_cast<uint32_t>(k[8]) << 24; [[fallthrough]];
559  case 8 : b += k[7]; [[fallthrough]];
560  case 7 : b += static_cast<uint32_t>(k[6]) << 8; [[fallthrough]];
561  case 6 : b += static_cast<uint32_t>(k[5]) << 16; [[fallthrough]];
562  case 5 : b += static_cast<uint32_t>(k[4]) << 24; [[fallthrough]];
563  case 4 : a += k[3]; [[fallthrough]];
564  case 3 : a += static_cast<uint32_t>(k[2]) << 8; [[fallthrough]];
565  case 2 : a += static_cast<uint32_t>(k[1]) << 16; [[fallthrough]];
566  case 1 : a += static_cast<uint32_t>(k[0]) << 24; [[fallthrough]];
567  break;
568  case 0 : return c;
569  }
570  }
571 
572  bitMixerFinal(a,b,c);
573  return c;
574 }
575 #endif
576 
577 
578 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
579 
580 unsigned Foam::Hasher
581 (
582  const void *key,
583  size_t length,
584  unsigned initval
585 )
586 {
587 #if defined (WM_BIG_ENDIAN)
588  return jenkins_hashbig(key, length, initval);
589 #elif defined (WM_LITTLE_ENDIAN)
590  return jenkins_hashlittle(key, length, initval);
591 #else
592  #error "Cannot determine WM_BIG_ENDIAN or WM_LITTLE_ENDIAN."
593 #endif
594 }
595 
596 
597 // ************************************************************************* //
HasherInt.H
Optimized hashing functions.
bitMixer
#define bitMixer(a, b, c)
Definition: Hasher.C:115
Foam::Hasher
unsigned Hasher(const void *data, size_t len, unsigned seed=0)
Bob Jenkins's 96-bit mixer hashing function (lookup3)
Definition: Hasher.C:581
Foam::HasherDual
unsigned HasherDual(const uint32_t *data, size_t length, unsigned &hash1, unsigned &hash2)
An optimized version of Hasher, returning dual hash values.
Foam::glTF::key
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
Definition: foamGltfBase.H:108
bitMixerFinal
#define bitMixerFinal(a, b, c)
Definition: Hasher.C:150
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Hasher.H
Miscellaneous hashing functions, mostly from Bob Jenkins.
Foam::HasherInt
unsigned HasherInt(const uint32_t *data, size_t length, unsigned seed=0)
An optimized version of Hasher.
Definition: Hasher.C:176
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
endian.H
Help with architecture-specific aspects.