xref: /petsc/src/sys/tests/ex63.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
18e89d99cSJDBetteridge 
28e89d99cSJDBetteridge static char help[] = "Tests `GarbageKeyAllReduceIntersect_Private()` in parallel\n\n";
38e89d99cSJDBetteridge 
48e89d99cSJDBetteridge #include <petscsys.h>
58e89d99cSJDBetteridge #include <petsc/private/garbagecollector.h>
68e89d99cSJDBetteridge 
78e89d99cSJDBetteridge /* This program tests `GarbageKeyAllReduceIntersect_Private()`.
88e89d99cSJDBetteridge    To test this routine in parallel, the sieve of Eratosthenes is
98e89d99cSJDBetteridge    implemented.
108e89d99cSJDBetteridge */
118e89d99cSJDBetteridge 
128e89d99cSJDBetteridge /* Populate an array with Prime numbers <= n.
138e89d99cSJDBetteridge    Primes are generated using trial division
148e89d99cSJDBetteridge */
158e89d99cSJDBetteridge PetscErrorCode Prime(PetscInt64 **set, PetscInt n)
168e89d99cSJDBetteridge {
178e89d99cSJDBetteridge   size_t      overestimate;
188e89d99cSJDBetteridge   PetscBool   is_prime;
198e89d99cSJDBetteridge   PetscInt64  ii, jj, count = 0;
208e89d99cSJDBetteridge   PetscInt64 *prime;
218e89d99cSJDBetteridge 
228e89d99cSJDBetteridge   PetscFunctionBeginUser;
238e89d99cSJDBetteridge   /* There will be fewer than ceil(1.26 * n/log(n)) primes <= n */
248e89d99cSJDBetteridge   overestimate = (size_t)PetscCeilReal(((PetscReal)n) * 1.26 / PetscLogReal((PetscReal)n));
258e89d99cSJDBetteridge   PetscCall(PetscMalloc1(overestimate, &prime));
268e89d99cSJDBetteridge   for (ii = 2; ii < n + 1; ii++) {
278e89d99cSJDBetteridge     is_prime = PETSC_TRUE;
288e89d99cSJDBetteridge     for (jj = 2; jj <= PetscFloorReal(PetscSqrtReal(ii)); jj++) {
298e89d99cSJDBetteridge       if (ii % jj == 0) {
308e89d99cSJDBetteridge         is_prime = PETSC_FALSE;
318e89d99cSJDBetteridge         break;
328e89d99cSJDBetteridge       }
338e89d99cSJDBetteridge     }
348e89d99cSJDBetteridge     if (is_prime) {
358e89d99cSJDBetteridge       prime[count] = ii;
368e89d99cSJDBetteridge       count++;
378e89d99cSJDBetteridge     }
388e89d99cSJDBetteridge   }
398e89d99cSJDBetteridge 
408e89d99cSJDBetteridge   PetscCall(PetscMalloc1((size_t)count + 1, set));
418e89d99cSJDBetteridge   (*set)[0] = count;
428e89d99cSJDBetteridge   for (ii = 1; ii < count + 1; ii++) { (*set)[ii] = prime[ii - 1]; }
438e89d99cSJDBetteridge   PetscCall(PetscFree(prime));
44*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
458e89d99cSJDBetteridge }
468e89d99cSJDBetteridge 
478e89d99cSJDBetteridge /* Print out the contents of a set */
488e89d99cSJDBetteridge PetscErrorCode PrintSet(MPI_Comm comm, PetscInt64 *set)
498e89d99cSJDBetteridge {
508e89d99cSJDBetteridge   char     text[64];
518e89d99cSJDBetteridge   PetscInt ii;
528e89d99cSJDBetteridge 
538e89d99cSJDBetteridge   PetscFunctionBeginUser;
548e89d99cSJDBetteridge   PetscCall(PetscSynchronizedPrintf(comm, "["));
558e89d99cSJDBetteridge   for (ii = 1; ii <= (PetscInt)set[0]; ii++) {
568e89d99cSJDBetteridge     PetscCall(PetscFormatConvert(" %" PetscInt64_FMT ",", text));
578e89d99cSJDBetteridge     PetscCall(PetscSynchronizedPrintf(comm, text, set[ii]));
588e89d99cSJDBetteridge   }
598e89d99cSJDBetteridge   PetscCall(PetscSynchronizedPrintf(comm, "]\n"));
608e89d99cSJDBetteridge   PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
61*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
628e89d99cSJDBetteridge }
638e89d99cSJDBetteridge 
648e89d99cSJDBetteridge /* Check set equality */
658e89d99cSJDBetteridge PetscErrorCode AssertSetsEqual(PetscInt64 *set, PetscInt64 *true_set)
668e89d99cSJDBetteridge {
678e89d99cSJDBetteridge   PetscInt ii;
688e89d99cSJDBetteridge 
698e89d99cSJDBetteridge   PetscFunctionBeginUser;
708e89d99cSJDBetteridge   PetscAssert((set[0] == true_set[0]), PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "Sets of different sizes");
718e89d99cSJDBetteridge   for (ii = 1; ii < set[0] + 1; ii++) { PetscAssert((set[ii] == true_set[ii]), PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "Sets are different"); }
72*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
738e89d99cSJDBetteridge }
748e89d99cSJDBetteridge 
758e89d99cSJDBetteridge /* Parallel implementation of the sieve of Eratosthenes */
768e89d99cSJDBetteridge PetscErrorCode test_sieve(MPI_Comm comm)
778e89d99cSJDBetteridge {
788e89d99cSJDBetteridge   PetscInt64  ii, local_p, maximum, n;
798e89d99cSJDBetteridge   PetscInt64 *local_set, *cursor, *bootstrap_primes, *truth;
808e89d99cSJDBetteridge   PetscMPIInt size, rank;
818e89d99cSJDBetteridge   PetscReal   x;
828e89d99cSJDBetteridge 
838e89d99cSJDBetteridge   PetscFunctionBeginUser;
848e89d99cSJDBetteridge   PetscCallMPI(MPI_Comm_rank(comm, &rank));
858e89d99cSJDBetteridge   PetscCallMPI(MPI_Comm_size(comm, &size));
868e89d99cSJDBetteridge 
878e89d99cSJDBetteridge   /* There should be at least `size + 1` primes smaller than
888e89d99cSJDBetteridge      (size + 1)*(log(size + 1) + log(log(size + 1))
898e89d99cSJDBetteridge     once `size` >=6
908e89d99cSJDBetteridge     This is sufficient for each rank to create its own sieve based off
918e89d99cSJDBetteridge     a different prime and calculate the size of the sieve.
928e89d99cSJDBetteridge   */
938e89d99cSJDBetteridge   x = (PetscReal)(size > 6) ? size + 1 : 7;
948e89d99cSJDBetteridge   x = x * (PetscLogReal(x) + PetscLogReal(PetscLogReal(x)));
958e89d99cSJDBetteridge   PetscCall(Prime(&bootstrap_primes, PetscCeilReal(x)));
968e89d99cSJDBetteridge 
978e89d99cSJDBetteridge   /* Calculate the maximum possible prime, select a prime number for
988e89d99cSJDBetteridge      each rank and allocate memorty for the sieve
998e89d99cSJDBetteridge   */
1008e89d99cSJDBetteridge   maximum = bootstrap_primes[size + 1] * bootstrap_primes[size + 1] - 1;
1018e89d99cSJDBetteridge   local_p = bootstrap_primes[rank + 1];
1028e89d99cSJDBetteridge   n       = maximum - local_p - (maximum / local_p) + 1 + rank + 1;
1038e89d99cSJDBetteridge   PetscCall(PetscMalloc1(n + 1, &local_set));
1048e89d99cSJDBetteridge 
1058e89d99cSJDBetteridge   /* Populate the sieve first with all primes <= `local_p`, followed by
1068e89d99cSJDBetteridge      all integers that are not a multiple of `local_p`
1078e89d99cSJDBetteridge   */
1088e89d99cSJDBetteridge   local_set[0] = n;
1098e89d99cSJDBetteridge   cursor       = &local_set[1];
1108e89d99cSJDBetteridge   for (ii = 0; ii < rank + 1; ii++) {
1118e89d99cSJDBetteridge     *cursor = bootstrap_primes[ii + 1];
1128e89d99cSJDBetteridge     cursor++;
1138e89d99cSJDBetteridge   }
1148e89d99cSJDBetteridge   for (ii = local_p + 1; ii <= maximum; ii++) {
1158e89d99cSJDBetteridge     if (ii % local_p != 0) {
1168e89d99cSJDBetteridge       *cursor = ii;
1178e89d99cSJDBetteridge       cursor++;
1188e89d99cSJDBetteridge     }
1198e89d99cSJDBetteridge   }
1208e89d99cSJDBetteridge   PetscCall(PetscPrintf(comm, "SIEVES:\n"));
1218e89d99cSJDBetteridge   PetscCall(PrintSet(comm, local_set));
1228e89d99cSJDBetteridge 
1238e89d99cSJDBetteridge   PetscCall(PetscFree(bootstrap_primes));
1248e89d99cSJDBetteridge 
1258e89d99cSJDBetteridge   /* Perform the intersection, testing parallel intersection routine */
1268e89d99cSJDBetteridge   PetscCall(GarbageKeyAllReduceIntersect_Private(PETSC_COMM_WORLD, &local_set[1], (PetscInt *)&local_set[0]));
1278e89d99cSJDBetteridge 
1288e89d99cSJDBetteridge   PetscCall(PetscPrintf(comm, "INTERSECTION:\n"));
1298e89d99cSJDBetteridge   PetscCall(PrintSet(comm, local_set));
1308e89d99cSJDBetteridge 
1318e89d99cSJDBetteridge   PetscCall(Prime(&truth, maximum));
1328e89d99cSJDBetteridge   PetscCall(PetscPrintf(comm, "TRUTH:\n"));
1338e89d99cSJDBetteridge   PetscCall(PrintSet(comm, truth));
1348e89d99cSJDBetteridge 
1358e89d99cSJDBetteridge   /* Assert the intersection corresponds to primes calculated using
1368e89d99cSJDBetteridge      trial division
1378e89d99cSJDBetteridge   */
1388e89d99cSJDBetteridge   PetscCall(AssertSetsEqual(local_set, truth));
1398e89d99cSJDBetteridge 
1408e89d99cSJDBetteridge   PetscCall(PetscFree(local_set));
1418e89d99cSJDBetteridge   PetscCall(PetscFree(truth));
142*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1438e89d99cSJDBetteridge }
1448e89d99cSJDBetteridge 
1458e89d99cSJDBetteridge /* Main executes the individual tests in a predefined order */
1468e89d99cSJDBetteridge int main(int argc, char **argv)
1478e89d99cSJDBetteridge {
1488e89d99cSJDBetteridge   PetscFunctionBeginUser;
1498e89d99cSJDBetteridge   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
1508e89d99cSJDBetteridge 
1518e89d99cSJDBetteridge   PetscCall(test_sieve(PETSC_COMM_WORLD));
1528e89d99cSJDBetteridge 
1538e89d99cSJDBetteridge   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "ALL PASSED\n"));
1548e89d99cSJDBetteridge   PetscCall(PetscFinalize());
1558e89d99cSJDBetteridge   return 0;
1568e89d99cSJDBetteridge }
1578e89d99cSJDBetteridge 
1588e89d99cSJDBetteridge /*TEST
1598e89d99cSJDBetteridge 
1608e89d99cSJDBetteridge    testset:
1618e89d99cSJDBetteridge      test:
1628e89d99cSJDBetteridge        nsize: 2
1638e89d99cSJDBetteridge        suffix: 2
1648e89d99cSJDBetteridge        output_file: output/ex63_2.out
1658e89d99cSJDBetteridge      test:
1668e89d99cSJDBetteridge        nsize: 3
1678e89d99cSJDBetteridge        suffix: 3
1688e89d99cSJDBetteridge        output_file: output/ex63_3.out
1698e89d99cSJDBetteridge      test:
1708e89d99cSJDBetteridge        nsize: 4
1718e89d99cSJDBetteridge        suffix: 4
1728e89d99cSJDBetteridge        output_file: output/ex63_4.out
1738e89d99cSJDBetteridge 
1748e89d99cSJDBetteridge TEST*/
175