xref: /petsc/src/sys/classes/random/tests/ex3.c (revision d0e6bf2ad94dcc89b258ce16c7987200a4714786)
1c4762a1bSJed Brown static char help[] = "Run Birthday Spacing Tests for PetscRandom.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscsys.h>
4c4762a1bSJed Brown #include <petscviewer.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown /* L'Ecuyer & Simard, 2001.
7c4762a1bSJed Brown  * "On the performance of birthday spacings tests with certain families of random number generators"
8c4762a1bSJed Brown  * https://doi.org/10.1016/S0378-4754(00)00253-6
9c4762a1bSJed Brown  */
10c4762a1bSJed Brown 
PoissonTailProbability(PetscReal lambda,PetscInt Y,PetscReal * prob)11d71ae5a4SJacob Faibussowitsch static PetscErrorCode PoissonTailProbability(PetscReal lambda, PetscInt Y, PetscReal *prob)
12d71ae5a4SJacob Faibussowitsch {
13c4762a1bSJed Brown   PetscReal p = 1.;
14c4762a1bSJed Brown   PetscReal logLambda;
15c4762a1bSJed Brown   PetscInt  i;
16c4762a1bSJed Brown   PetscReal logFact = 0.;
17c4762a1bSJed Brown 
18c4762a1bSJed Brown   PetscFunctionBegin;
19c4762a1bSJed Brown   logLambda = PetscLogReal(lambda);
20c4762a1bSJed Brown   for (i = 0; i < Y; i++) {
21c4762a1bSJed Brown     PetscReal exponent = -lambda + i * logLambda - logFact;
22c4762a1bSJed Brown 
23c4762a1bSJed Brown     logFact += PetscLogReal(i + 1);
24c4762a1bSJed Brown 
25c4762a1bSJed Brown     p -= PetscExpReal(exponent);
26c4762a1bSJed Brown   }
27c4762a1bSJed Brown   *prob = p;
283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29c4762a1bSJed Brown }
30c4762a1bSJed Brown 
main(int argc,char ** argv)31d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
32d71ae5a4SJacob Faibussowitsch {
33c4762a1bSJed Brown   PetscMPIInt size;
34c4762a1bSJed Brown   PetscInt    log2d, log2n, t, N, Y;
35c4762a1bSJed Brown   PetscInt64  d, k, *X;
366497c311SBarry Smith   PetscInt    n, i;
37c4762a1bSJed Brown   PetscReal   lambda, p;
38c4762a1bSJed Brown   PetscRandom random;
39c4762a1bSJed Brown 
40327415f7SBarry Smith   PetscFunctionBeginUser;
419566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
42c4762a1bSJed Brown   t     = 8;
43c4762a1bSJed Brown   log2d = 7;
44c4762a1bSJed Brown   log2n = 20;
45d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Birthday spacing test parameters", "PetscRandom");
469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-t", "t, the dimension of the sample space", "ex3.c", t, &t, NULL));
479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-log2d", "The log of d, the number of bins per direction", "ex3.c", log2d, &log2d, NULL));
489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-log2n", "The log of n, the number of samples per process", "ex3.c", log2n, &log2n, NULL));
49d0609cedSBarry Smith   PetscOptionsEnd();
50c4762a1bSJed Brown 
51cc73adaaSBarry Smith   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);
52c4762a1bSJed Brown   d = ((PetscInt64)1) << log2d;
53c4762a1bSJed Brown   k = ((PetscInt64)1) << (log2d * t);
54cc73adaaSBarry Smith   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);
556497c311SBarry Smith   n = (PetscInt)(((size_t)1) << log2n);
569566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
57c4762a1bSJed Brown   N      = size;
58*57508eceSPierre Jolivet   lambda = PetscPowRealInt(2., 3 * log2n - (2 + log2d * t));
59c4762a1bSJed Brown 
606497c311SBarry Smith   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Generating %" PetscInt_FMT " 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));
619566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Expected spacing collisions per process %g (%g total).\n", (double)lambda, (double)(N * lambda)));
629566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &random));
639566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(random));
649566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(random, 0.0, 1.0));
659566063dSJacob Faibussowitsch   PetscCall(PetscRandomView(random, PETSC_VIEWER_STDOUT_WORLD));
669566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &X));
67c4762a1bSJed Brown   for (i = 0; i < n; i++) {
68c4762a1bSJed Brown     PetscInt   j;
69c4762a1bSJed Brown     PetscInt64 bin  = 0;
70c4762a1bSJed Brown     PetscInt64 mult = 1;
71c4762a1bSJed Brown 
72c4762a1bSJed Brown     for (j = 0; j < t; j++, mult *= d) {
73c4762a1bSJed Brown       PetscReal  x;
74c4762a1bSJed Brown       PetscInt64 slot;
75c4762a1bSJed Brown 
769566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValueReal(random, &x));
77c4762a1bSJed Brown       /* worried about precision here */
78c4762a1bSJed Brown       slot = (PetscInt64)(x * d);
79c4762a1bSJed Brown       bin += mult * slot;
80c4762a1bSJed Brown     }
8108401ef6SPierre Jolivet     PetscCheck(bin < k, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Generated point in bin %" PetscInt64_FMT ", but only %" PetscInt64_FMT " bins", bin, k);
82c4762a1bSJed Brown     X[i] = bin;
83c4762a1bSJed Brown   }
84c4762a1bSJed Brown 
857a43aa34SPierre Jolivet   PetscCall(PetscSortInt64(n, X));
86ad540459SPierre Jolivet   for (i = 0; i < n - 1; i++) X[i] = X[i + 1] - X[i];
877a43aa34SPierre Jolivet   PetscCall(PetscSortInt64(n - 1, X));
88ad540459SPierre Jolivet   for (i = 0, Y = 0; i < n - 2; i++) Y += (X[i + 1] == X[i]);
89c4762a1bSJed Brown 
90462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &Y, 1, MPIU_INT, MPI_SUM, MPI_COMM_WORLD));
919566063dSJacob Faibussowitsch   PetscCall(PoissonTailProbability(N * lambda, Y, &p));
92da81f932SPierre Jolivet   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT " total collisions counted: that many or more should occur with probability %g.\n", Y, (double)p));
93c4762a1bSJed Brown 
949566063dSJacob Faibussowitsch   PetscCall(PetscFree(X));
959566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&random));
969566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
97b122ec5aSJacob Faibussowitsch   return 0;
98c4762a1bSJed Brown }
99c4762a1bSJed Brown 
100c4762a1bSJed Brown /*TEST
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   test:
103c4762a1bSJed Brown     args: -t 4 -log2d 7 -log2n 10
104c4762a1bSJed Brown 
105c4762a1bSJed Brown TEST*/
106