xref: /petsc/src/vec/is/sf/tests/ex17.c (revision 061e922f3926be00487707c73b78dd3d40309129)
181b34b04SJunchao Zhang static const char help[] = "Test PetscSF with MPI large count (more than 2 billion elements in messages)\n\n";
281b34b04SJunchao Zhang 
381b34b04SJunchao Zhang #include <petscsys.h>
481b34b04SJunchao Zhang #include <petscsf.h>
581b34b04SJunchao Zhang 
main(int argc,char ** argv)6*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
7*d71ae5a4SJacob Faibussowitsch {
881b34b04SJunchao Zhang   PetscSF            sf;
981b34b04SJunchao Zhang   PetscInt           i, nroots, nleaves;
1081b34b04SJunchao Zhang   PetscInt           n       = (1ULL << 31) + 1024; /* a little over 2G elements */
1181b34b04SJunchao Zhang   PetscSFNode       *iremote = NULL;
1281b34b04SJunchao Zhang   PetscMPIInt        rank, size;
1381b34b04SJunchao Zhang   char              *rootdata = NULL, *leafdata = NULL;
1481b34b04SJunchao Zhang   Vec                x, y;
1581b34b04SJunchao Zhang   VecScatter         vscat;
1681b34b04SJunchao Zhang   PetscInt           rstart, rend;
1781b34b04SJunchao Zhang   IS                 ix;
1881b34b04SJunchao Zhang   const PetscScalar *xv;
1981b34b04SJunchao Zhang 
20327415f7SBarry Smith   PetscFunctionBeginUser;
219566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
2408401ef6SPierre Jolivet   PetscCheck(size == 2, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "The test can only run with two MPI ranks");
2581b34b04SJunchao Zhang 
2681b34b04SJunchao Zhang   /* Test PetscSF */
279566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sf));
289566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sf));
2981b34b04SJunchao Zhang 
30c5853193SPierre Jolivet   if (rank == 0) {
3181b34b04SJunchao Zhang     nroots  = n;
3281b34b04SJunchao Zhang     nleaves = 0;
3381b34b04SJunchao Zhang   } else {
3481b34b04SJunchao Zhang     nroots  = 0;
3581b34b04SJunchao Zhang     nleaves = n;
369566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves, &iremote));
3781b34b04SJunchao Zhang     for (i = 0; i < nleaves; i++) {
3881b34b04SJunchao Zhang       iremote[i].rank  = 0;
3981b34b04SJunchao Zhang       iremote[i].index = i;
4081b34b04SJunchao Zhang     }
4181b34b04SJunchao Zhang   }
429566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(sf, nroots, nleaves, NULL, PETSC_COPY_VALUES, iremote, PETSC_COPY_VALUES));
439566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nroots, &rootdata, nleaves, &leafdata));
44c5853193SPierre Jolivet   if (rank == 0) {
4581b34b04SJunchao Zhang     memset(rootdata, 11, nroots);
4681b34b04SJunchao Zhang     rootdata[nroots - 1] = 12; /* Use a different value at the end */
4781b34b04SJunchao Zhang   }
4881b34b04SJunchao Zhang 
499566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sf, MPI_SIGNED_CHAR, rootdata, leafdata, MPI_REPLACE)); /* rank 0->1, bcast rootdata to leafdata */
509566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sf, MPI_SIGNED_CHAR, rootdata, leafdata, MPI_REPLACE));
519566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(sf, MPI_SIGNED_CHAR, leafdata, rootdata, MPI_SUM)); /* rank 1->0, add leafdata to rootdata */
529566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd(sf, MPI_SIGNED_CHAR, leafdata, rootdata, MPI_SUM));
53c5853193SPierre Jolivet   PetscCheck(rank != 0 || (rootdata[0] == 22 && rootdata[nroots - 1] == 24), PETSC_COMM_SELF, PETSC_ERR_PLIB, "SF: wrong results");
5481b34b04SJunchao Zhang 
559566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rootdata, leafdata));
569566063dSJacob Faibussowitsch   PetscCall(PetscFree(iremote));
579566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sf));
5881b34b04SJunchao Zhang 
5981b34b04SJunchao Zhang   /* Test VecScatter */
609566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
619566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &y));
629566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, rank == 0 ? n : 64, PETSC_DECIDE));
639566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(y, rank == 0 ? 64 : n, PETSC_DECIDE));
649566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
659566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(y));
6681b34b04SJunchao Zhang 
679566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(x, &rstart, &rend));
689566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, rend - rstart, rstart, 1, &ix));
699566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x, ix, y, ix, &vscat));
7081b34b04SJunchao Zhang 
719566063dSJacob Faibussowitsch   PetscCall(VecSet(x, 3.0));
729566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
739566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(vscat, x, y, INSERT_VALUES, SCATTER_FORWARD));
7481b34b04SJunchao Zhang 
759566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(vscat, y, x, ADD_VALUES, SCATTER_REVERSE));
769566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(vscat, y, x, ADD_VALUES, SCATTER_REVERSE));
7781b34b04SJunchao Zhang 
789566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xv));
7908401ef6SPierre Jolivet   PetscCheck(xv[0] == 6.0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "VecScatter: wrong results");
809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xv));
819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
829566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
839566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&vscat));
849566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ix));
8581b34b04SJunchao Zhang 
869566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
87b122ec5aSJacob Faibussowitsch   return 0;
8881b34b04SJunchao Zhang }
8981b34b04SJunchao Zhang 
9081b34b04SJunchao Zhang /**TEST
9181b34b04SJunchao Zhang    test:
9281b34b04SJunchao Zhang      requires: defined(PETSC_HAVE_MPI_LARGE_COUNT) defined(PETSC_USE_64BIT_INDICES)
9381b34b04SJunchao Zhang      TODO: need a machine with big memory (~150GB) to run the test
9481b34b04SJunchao Zhang      nsize: 2
9581b34b04SJunchao Zhang      args: -sf_type {{basic neighbor}}
9681b34b04SJunchao Zhang 
9781b34b04SJunchao Zhang TEST**/
98