xref: /petsc/src/sys/classes/random/tests/ex3.c (revision f6dfbefd03961ab3be6f06be75c96cbf27a49667)
1 
2 static char help[] = "Run Birthday Spacing Tests for PetscRandom.\n\n";
3 
4 #include <petscsys.h>
5 #include <petscviewer.h>
6 
7 /* L'Ecuyer & Simard, 2001.
8  * "On the performance of birthday spacings tests with certain families of random number generators"
9  * https://doi.org/10.1016/S0378-4754(00)00253-6
10  */
11 
12 static int PetscInt64Compare(const void *a, const void *b) {
13   PetscInt64 A = *((const PetscInt64 *)a);
14   PetscInt64 B = *((const PetscInt64 *)b);
15   if (A < B) return -1;
16   if (A == B) return 0;
17   return 1;
18 }
19 
20 static PetscErrorCode PoissonTailProbability(PetscReal lambda, PetscInt Y, PetscReal *prob) {
21   PetscReal p = 1.;
22   PetscReal logLambda;
23   PetscInt  i;
24   PetscReal logFact = 0.;
25 
26   PetscFunctionBegin;
27   logLambda = PetscLogReal(lambda);
28   for (i = 0; i < Y; i++) {
29     PetscReal exponent = -lambda + i * logLambda - logFact;
30 
31     logFact += PetscLogReal(i + 1);
32 
33     p -= PetscExpReal(exponent);
34   }
35   *prob = p;
36   PetscFunctionReturn(0);
37 }
38 
39 int main(int argc, char **argv) {
40   PetscMPIInt size;
41   PetscInt    log2d, log2n, t, N, Y;
42   PetscInt64  d, k, *X;
43   size_t      n, i;
44   PetscReal   lambda, p;
45   PetscRandom random;
46 
47   PetscFunctionBeginUser;
48   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
49   t     = 8;
50   log2d = 7;
51   log2n = 20;
52   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Birthday spacing test parameters", "PetscRandom");
53   PetscCall(PetscOptionsInt("-t", "t, the dimension of the sample space", "ex3.c", t, &t, NULL));
54   PetscCall(PetscOptionsInt("-log2d", "The log of d, the number of bins per direction", "ex3.c", log2d, &log2d, NULL));
55   PetscCall(PetscOptionsInt("-log2n", "The log of n, the number of samples per process", "ex3.c", log2n, &log2n, NULL));
56   PetscOptionsEnd();
57 
58   PetscCheck((size_t)log2d * t <= sizeof(PetscInt64) * 8 - 2, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "The number of bins (2^%" PetscInt_FMT ") is too big for PetscInt64.", log2d * t);
59   d = ((PetscInt64)1) << log2d;
60   k = ((PetscInt64)1) << (log2d * t);
61   PetscCheck((size_t)log2n <= sizeof(size_t) * 8 - 1, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "The number of samples per process (2^%" PetscInt_FMT ") is too big for size_t.", log2n);
62   n = ((size_t)1) << log2n;
63   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
64   N      = size;
65   lambda = PetscPowRealInt(2., (3 * log2n - (2 + log2d * t)));
66 
67   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Generating %zu samples (%g GB) per process in a %" PetscInt_FMT " dimensional space with %" PetscInt64_FMT " bins per dimension.\n", n, (n * 1.e-9) * sizeof(PetscInt64), t, d));
68   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Expected spacing collisions per process %g (%g total).\n", (double)lambda, (double)(N * lambda)));
69   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &random));
70   PetscCall(PetscRandomSetFromOptions(random));
71   PetscCall(PetscRandomSetInterval(random, 0.0, 1.0));
72   PetscCall(PetscRandomView(random, PETSC_VIEWER_STDOUT_WORLD));
73   PetscCall(PetscMalloc1(n, &X));
74   for (i = 0; i < n; i++) {
75     PetscInt   j;
76     PetscInt64 bin  = 0;
77     PetscInt64 mult = 1;
78 
79     for (j = 0; j < t; j++, mult *= d) {
80       PetscReal  x;
81       PetscInt64 slot;
82 
83       PetscCall(PetscRandomGetValueReal(random, &x));
84       /* worried about precision here */
85       slot = (PetscInt64)(x * d);
86       bin += mult * slot;
87     }
88     PetscCheck(bin < k, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Generated point in bin %" PetscInt64_FMT ", but only %" PetscInt64_FMT " bins", bin, k);
89     X[i] = bin;
90   }
91 
92   qsort(X, n, sizeof(PetscInt64), PetscInt64Compare);
93   for (i = 0; i < n - 1; i++) X[i] = X[i + 1] - X[i];
94   qsort(X, n - 1, sizeof(PetscInt64), PetscInt64Compare);
95   for (i = 0, Y = 0; i < n - 2; i++) Y += (X[i + 1] == X[i]);
96 
97   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &Y, 1, MPIU_INT, MPI_SUM, MPI_COMM_WORLD));
98   PetscCall(PoissonTailProbability(N * lambda, Y, &p));
99   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT " total collisions counted: that many or more should occur with probabilty %g.\n", Y, (double)p));
100 
101   PetscCall(PetscFree(X));
102   PetscCall(PetscRandomDestroy(&random));
103   PetscCall(PetscFinalize());
104   return 0;
105 }
106 
107 /*TEST
108 
109   test:
110     args: -t 4 -log2d 7 -log2n 10
111 
112 TEST*/
113