xref: /petsc/src/vec/is/sf/tests/ex17.c (revision 061e922f3926be00487707c73b78dd3d40309129) !
1 static const char help[] = "Test PetscSF with MPI large count (more than 2 billion elements in messages)\n\n";
2 
3 #include <petscsys.h>
4 #include <petscsf.h>
5 
main(int argc,char ** argv)6 int main(int argc, char **argv)
7 {
8   PetscSF            sf;
9   PetscInt           i, nroots, nleaves;
10   PetscInt           n       = (1ULL << 31) + 1024; /* a little over 2G elements */
11   PetscSFNode       *iremote = NULL;
12   PetscMPIInt        rank, size;
13   char              *rootdata = NULL, *leafdata = NULL;
14   Vec                x, y;
15   VecScatter         vscat;
16   PetscInt           rstart, rend;
17   IS                 ix;
18   const PetscScalar *xv;
19 
20   PetscFunctionBeginUser;
21   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
22   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
23   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
24   PetscCheck(size == 2, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "The test can only run with two MPI ranks");
25 
26   /* Test PetscSF */
27   PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sf));
28   PetscCall(PetscSFSetFromOptions(sf));
29 
30   if (rank == 0) {
31     nroots  = n;
32     nleaves = 0;
33   } else {
34     nroots  = 0;
35     nleaves = n;
36     PetscCall(PetscMalloc1(nleaves, &iremote));
37     for (i = 0; i < nleaves; i++) {
38       iremote[i].rank  = 0;
39       iremote[i].index = i;
40     }
41   }
42   PetscCall(PetscSFSetGraph(sf, nroots, nleaves, NULL, PETSC_COPY_VALUES, iremote, PETSC_COPY_VALUES));
43   PetscCall(PetscMalloc2(nroots, &rootdata, nleaves, &leafdata));
44   if (rank == 0) {
45     memset(rootdata, 11, nroots);
46     rootdata[nroots - 1] = 12; /* Use a different value at the end */
47   }
48 
49   PetscCall(PetscSFBcastBegin(sf, MPI_SIGNED_CHAR, rootdata, leafdata, MPI_REPLACE)); /* rank 0->1, bcast rootdata to leafdata */
50   PetscCall(PetscSFBcastEnd(sf, MPI_SIGNED_CHAR, rootdata, leafdata, MPI_REPLACE));
51   PetscCall(PetscSFReduceBegin(sf, MPI_SIGNED_CHAR, leafdata, rootdata, MPI_SUM)); /* rank 1->0, add leafdata to rootdata */
52   PetscCall(PetscSFReduceEnd(sf, MPI_SIGNED_CHAR, leafdata, rootdata, MPI_SUM));
53   PetscCheck(rank != 0 || (rootdata[0] == 22 && rootdata[nroots - 1] == 24), PETSC_COMM_SELF, PETSC_ERR_PLIB, "SF: wrong results");
54 
55   PetscCall(PetscFree2(rootdata, leafdata));
56   PetscCall(PetscFree(iremote));
57   PetscCall(PetscSFDestroy(&sf));
58 
59   /* Test VecScatter */
60   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
61   PetscCall(VecCreate(PETSC_COMM_WORLD, &y));
62   PetscCall(VecSetSizes(x, rank == 0 ? n : 64, PETSC_DECIDE));
63   PetscCall(VecSetSizes(y, rank == 0 ? 64 : n, PETSC_DECIDE));
64   PetscCall(VecSetFromOptions(x));
65   PetscCall(VecSetFromOptions(y));
66 
67   PetscCall(VecGetOwnershipRange(x, &rstart, &rend));
68   PetscCall(ISCreateStride(PETSC_COMM_SELF, rend - rstart, rstart, 1, &ix));
69   PetscCall(VecScatterCreate(x, ix, y, ix, &vscat));
70 
71   PetscCall(VecSet(x, 3.0));
72   PetscCall(VecScatterBegin(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
73   PetscCall(VecScatterEnd(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
74 
75   PetscCall(VecScatterBegin(vscat, y, x, ADD_VALUES, SCATTER_REVERSE));
76   PetscCall(VecScatterEnd(vscat, y, x, ADD_VALUES, SCATTER_REVERSE));
77 
78   PetscCall(VecGetArrayRead(x, &xv));
79   PetscCheck(xv[0] == 6.0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "VecScatter: wrong results");
80   PetscCall(VecRestoreArrayRead(x, &xv));
81   PetscCall(VecDestroy(&x));
82   PetscCall(VecDestroy(&y));
83   PetscCall(VecScatterDestroy(&vscat));
84   PetscCall(ISDestroy(&ix));
85 
86   PetscCall(PetscFinalize());
87   return 0;
88 }
89 
90 /**TEST
91    test:
92      requires: defined(PETSC_HAVE_MPI_LARGE_COUNT) defined(PETSC_USE_64BIT_INDICES)
93      TODO: need a machine with big memory (~150GB) to run the test
94      nsize: 2
95      args: -sf_type {{basic neighbor}}
96 
97 TEST**/
98