sobol_gold.c 5.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2010-2011 Université de Bordeaux 1
  4. *
  5. * StarPU is free software; you can redistribute it and/or modify
  6. * it under the terms of the GNU Lesser General Public License as published by
  7. * the Free Software Foundation; either version 2.1 of the License, or (at
  8. * your option) any later version.
  9. *
  10. * StarPU is distributed in the hope that it will be useful, but
  11. * WITHOUT ANY WARRANTY; without even the implied warranty of
  12. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
  13. *
  14. * See the GNU Lesser General Public License in COPYING.LGPL for more details.
  15. */
  16. /*
  17. * Copyright 1993-2009 NVIDIA Corporation. All rights reserved.
  18. *
  19. * NVIDIA Corporation and its licensors retain all intellectual property and
  20. * proprietary rights in and to this software and related documentation and
  21. * any modifications thereto. Any use, reproduction, disclosure, or distribution
  22. * of this software and related documentation without an express license
  23. * agreement from NVIDIA Corporation is strictly prohibited.
  24. *
  25. */
  26. /*
  27. * Portions Copyright (c) 1993-2009 NVIDIA Corporation. All rights reserved.
  28. * Portions Copyright (c) 2009 Mike Giles, Oxford University. All rights reserved.
  29. * Portions Copyright (c) 2008 Frances Y. Kuo and Stephen Joe. All rights reserved.
  30. *
  31. * Sobol Quasi-random Number Generator example
  32. *
  33. * Based on CUDA code submitted by Mike Giles, Oxford University, United Kingdom
  34. * http://people.maths.ox.ac.uk/~gilesm/
  35. *
  36. * and C code developed by Stephen Joe, University of Waikato, New Zealand
  37. * and Frances Kuo, University of New South Wales, Australia
  38. * http://web.maths.unsw.edu.au/~fkuo/sobol/
  39. *
  40. * For theoretical background see:
  41. *
  42. * P. Bratley and B.L. Fox.
  43. * Implementing Sobol's quasirandom sequence generator
  44. * http://portal.acm.org/citation.cfm?id=42288
  45. * ACM Trans. on Math. Software, 14(1):88-100, 1988
  46. *
  47. * S. Joe and F. Kuo.
  48. * Remark on algorithm 659: implementing Sobol's quasirandom sequence generator.
  49. * http://portal.acm.org/citation.cfm?id=641879
  50. * ACM Trans. on Math. Software, 29(1):49-57, 2003
  51. */
  52. #include <stdio.h>
  53. #include <stdlib.h>
  54. #include <math.h>
  55. #include <string.h>
  56. #include <strings.h>
  57. #include "sobol.h"
  58. #include "sobol_gold.h"
  59. #include "sobol_primitives.h"
  60. #define k_2powneg32 2.3283064E-10F
  61. /* Create the direction numbers, based on the primitive polynomials. */
  62. void initSobolDirectionVectors(int n_dimensions, unsigned int *directions)
  63. {
  64. unsigned int *v = directions;
  65. int dim;
  66. for (dim = 0 ; dim < n_dimensions ; dim++)
  67. {
  68. /* First dimension is a special case */
  69. if (dim == 0)
  70. {
  71. int i;
  72. for (i = 0 ; i < n_directions ; i++)
  73. {
  74. /* All m's are 1 */
  75. v[i] = 1 << (31 - i);
  76. }
  77. }
  78. else
  79. {
  80. int d = sobol_primitives[dim].degree;
  81. /* The first direction numbers (up to the degree of the polynomial)
  82. are simply v[i] = m[i] / 2^i (stored in Q0.32 format) */
  83. int i;
  84. for (i = 0 ; i < d ; i++)
  85. {
  86. v[i] = sobol_primitives[dim].m[i] << (31 - i);
  87. }
  88. /* The remaining direction numbers are computed as described in
  89. the Bratley and Fox paper. */
  90. /* v[i] = a[1]v[i-1] ^ a[2]v[i-2] ^ ... ^ a[v-1]v[i-d+1] ^ v[i-d] ^ v[i-d]/2^d */
  91. for (i = d ; i < n_directions ; i++)
  92. {
  93. /* First do the v[i-d] ^ v[i-d]/2^d part */
  94. v[i] = v[i - d] ^ (v[i - d] >> d);
  95. /* Now do the a[1]v[i-1] ^ a[2]v[i-2] ^ ... part
  96. Note that the coefficients a[] are zero or one and for compactness in
  97. the input tables they are stored as bits of a single integer. To extract
  98. the relevant bit we use right shift and mask with 1.
  99. For example, for a 10 degree polynomial there are ten useful bits in a,
  100. so to get a[2] we need to right shift 7 times (to get the 8th bit into
  101. the LSB) and then mask with 1. */
  102. int j;
  103. for (j = 1 ; j < d ; j++)
  104. {
  105. v[i] ^= (((sobol_primitives[dim].a >> (d - 1 - j)) & 1) * v[i - j]);
  106. }
  107. }
  108. }
  109. v += n_directions;
  110. }
  111. }
  112. /* Reference model for generating Sobol numbers on the host */
  113. void sobolCPU(int n_vectors, int n_dimensions, unsigned int *directions, float *output)
  114. {
  115. unsigned int *v = directions;
  116. int d;
  117. for (d = 0 ; d < n_dimensions ; d++)
  118. {
  119. unsigned int X = 0;
  120. /* x[0] is zero (in all dimensions) */
  121. output[n_vectors * d] = 0.0;
  122. int i;
  123. for (i = 1 ; i < n_vectors ; i++)
  124. {
  125. /* x[i] = x[i-1] ^ v[c]
  126. where c is the index of the rightmost zero bit in i
  127. minus 1 (since C arrays count from zero)
  128. In the Bratley and Fox paper this is equation (**) */
  129. X ^= v[ffs(~(i - 1)) - 1];
  130. output[i + n_vectors * d] = (float)X * k_2powneg32;
  131. }
  132. v += n_directions;
  133. }
  134. }