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