1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Run Birthday Spacing Tests for PetscRandom.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscsys.h> 5c4762a1bSJed Brown #include <petscviewer.h> 6c4762a1bSJed Brown 7c4762a1bSJed Brown /* L'Ecuyer & Simard, 2001. 8c4762a1bSJed Brown * "On the performance of birthday spacings tests with certain families of random number generators" 9c4762a1bSJed Brown * https://doi.org/10.1016/S0378-4754(00)00253-6 10c4762a1bSJed Brown */ 11c4762a1bSJed Brown 129371c9d4SSatish Balay static int PetscInt64Compare(const void *a, const void *b) { 13c4762a1bSJed Brown PetscInt64 A = *((const PetscInt64 *)a); 14c4762a1bSJed Brown PetscInt64 B = *((const PetscInt64 *)b); 15c4762a1bSJed Brown if (A < B) return -1; 16c4762a1bSJed Brown if (A == B) return 0; 17c4762a1bSJed Brown return 1; 18c4762a1bSJed Brown } 19c4762a1bSJed Brown 209371c9d4SSatish Balay static PetscErrorCode PoissonTailProbability(PetscReal lambda, PetscInt Y, PetscReal *prob) { 21c4762a1bSJed Brown PetscReal p = 1.; 22c4762a1bSJed Brown PetscReal logLambda; 23c4762a1bSJed Brown PetscInt i; 24c4762a1bSJed Brown PetscReal logFact = 0.; 25c4762a1bSJed Brown 26c4762a1bSJed Brown PetscFunctionBegin; 27c4762a1bSJed Brown logLambda = PetscLogReal(lambda); 28c4762a1bSJed Brown for (i = 0; i < Y; i++) { 29c4762a1bSJed Brown PetscReal exponent = -lambda + i * logLambda - logFact; 30c4762a1bSJed Brown 31c4762a1bSJed Brown logFact += PetscLogReal(i + 1); 32c4762a1bSJed Brown 33c4762a1bSJed Brown p -= PetscExpReal(exponent); 34c4762a1bSJed Brown } 35c4762a1bSJed Brown *prob = p; 36c4762a1bSJed Brown PetscFunctionReturn(0); 37c4762a1bSJed Brown } 38c4762a1bSJed Brown 399371c9d4SSatish Balay int main(int argc, char **argv) { 40c4762a1bSJed Brown PetscMPIInt size; 41c4762a1bSJed Brown PetscInt log2d, log2n, t, N, Y; 42c4762a1bSJed Brown PetscInt64 d, k, *X; 43c4762a1bSJed Brown size_t n, i; 44c4762a1bSJed Brown PetscReal lambda, p; 45c4762a1bSJed Brown PetscRandom random; 46c4762a1bSJed Brown 47327415f7SBarry Smith PetscFunctionBeginUser; 489566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 49c4762a1bSJed Brown t = 8; 50c4762a1bSJed Brown log2d = 7; 51c4762a1bSJed Brown log2n = 20; 52d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Birthday spacing test parameters", "PetscRandom"); 539566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-t", "t, the dimension of the sample space", "ex3.c", t, &t, NULL)); 549566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-log2d", "The log of d, the number of bins per direction", "ex3.c", log2d, &log2d, NULL)); 559566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-log2n", "The log of n, the number of samples per process", "ex3.c", log2n, &log2n, NULL)); 56d0609cedSBarry Smith PetscOptionsEnd(); 57c4762a1bSJed Brown 58cc73adaaSBarry 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); 59c4762a1bSJed Brown d = ((PetscInt64)1) << log2d; 60c4762a1bSJed Brown k = ((PetscInt64)1) << (log2d * t); 61cc73adaaSBarry 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); 62c4762a1bSJed Brown n = ((size_t)1) << log2n; 639566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 64c4762a1bSJed Brown N = size; 65c4762a1bSJed Brown lambda = PetscPowRealInt(2., (3 * log2n - (2 + log2d * t))); 66c4762a1bSJed Brown 679566063dSJacob Faibussowitsch 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)); 689566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Expected spacing collisions per process %g (%g total).\n", (double)lambda, (double)(N * lambda))); 699566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &random)); 709566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(random)); 719566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(random, 0.0, 1.0)); 729566063dSJacob Faibussowitsch PetscCall(PetscRandomView(random, PETSC_VIEWER_STDOUT_WORLD)); 739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &X)); 74c4762a1bSJed Brown for (i = 0; i < n; i++) { 75c4762a1bSJed Brown PetscInt j; 76c4762a1bSJed Brown PetscInt64 bin = 0; 77c4762a1bSJed Brown PetscInt64 mult = 1; 78c4762a1bSJed Brown 79c4762a1bSJed Brown for (j = 0; j < t; j++, mult *= d) { 80c4762a1bSJed Brown PetscReal x; 81c4762a1bSJed Brown PetscInt64 slot; 82c4762a1bSJed Brown 839566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(random, &x)); 84c4762a1bSJed Brown /* worried about precision here */ 85c4762a1bSJed Brown slot = (PetscInt64)(x * d); 86c4762a1bSJed Brown bin += mult * slot; 87c4762a1bSJed Brown } 8808401ef6SPierre Jolivet PetscCheck(bin < k, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Generated point in bin %" PetscInt64_FMT ", but only %" PetscInt64_FMT " bins", bin, k); 89c4762a1bSJed Brown X[i] = bin; 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 92c4762a1bSJed Brown qsort(X, n, sizeof(PetscInt64), PetscInt64Compare); 93*ad540459SPierre Jolivet for (i = 0; i < n - 1; i++) X[i] = X[i + 1] - X[i]; 94c4762a1bSJed Brown qsort(X, n - 1, sizeof(PetscInt64), PetscInt64Compare); 95*ad540459SPierre Jolivet for (i = 0, Y = 0; i < n - 2; i++) Y += (X[i + 1] == X[i]); 96c4762a1bSJed Brown 979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &Y, 1, MPIU_INT, MPI_SUM, MPI_COMM_WORLD)); 989566063dSJacob Faibussowitsch PetscCall(PoissonTailProbability(N * lambda, Y, &p)); 999566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT " total collisions counted: that many or more should occur with probabilty %g.\n", Y, (double)p)); 100c4762a1bSJed Brown 1019566063dSJacob Faibussowitsch PetscCall(PetscFree(X)); 1029566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&random)); 1039566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 104b122ec5aSJacob Faibussowitsch return 0; 105c4762a1bSJed Brown } 106c4762a1bSJed Brown 107c4762a1bSJed Brown /*TEST 108c4762a1bSJed Brown 109c4762a1bSJed Brown test: 110c4762a1bSJed Brown args: -t 4 -log2d 7 -log2n 10 111c4762a1bSJed Brown 112c4762a1bSJed Brown TEST*/ 113