xref: /petsc/src/vec/is/sf/tests/ex17.c (revision 9566063d113dddea24716c546802770db7481bc0)
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 
681b34b04SJunchao Zhang int main(int argc,char **argv)
781b34b04SJunchao Zhang {
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 
20*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
21*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
22*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
232c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"The test can only run with two MPI ranks");
2481b34b04SJunchao Zhang 
2581b34b04SJunchao Zhang   /* Test PetscSF */
26*9566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PETSC_COMM_WORLD,&sf));
27*9566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sf));
2881b34b04SJunchao Zhang 
2981b34b04SJunchao Zhang   if (!rank) {
3081b34b04SJunchao Zhang     nroots  = n;
3181b34b04SJunchao Zhang     nleaves = 0;
3281b34b04SJunchao Zhang   } else {
3381b34b04SJunchao Zhang     nroots  = 0;
3481b34b04SJunchao Zhang     nleaves = n;
35*9566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves,&iremote));
3681b34b04SJunchao Zhang     for (i=0; i<nleaves; i++) {
3781b34b04SJunchao Zhang       iremote[i].rank  = 0;
3881b34b04SJunchao Zhang       iremote[i].index = i;
3981b34b04SJunchao Zhang     }
4081b34b04SJunchao Zhang   }
41*9566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(sf,nroots,nleaves,NULL,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES));
42*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nroots,&rootdata,nleaves,&leafdata));
4381b34b04SJunchao Zhang   if (!rank) {
4481b34b04SJunchao Zhang     memset(rootdata,11,nroots);
4581b34b04SJunchao Zhang     rootdata[nroots-1] = 12; /* Use a different value at the end */
4681b34b04SJunchao Zhang   }
4781b34b04SJunchao Zhang 
48*9566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sf,MPI_SIGNED_CHAR,rootdata,leafdata,MPI_REPLACE)); /* rank 0->1, bcast rootdata to leafdata */
49*9566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sf,MPI_SIGNED_CHAR,rootdata,leafdata,MPI_REPLACE));
50*9566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(sf,MPI_SIGNED_CHAR,leafdata,rootdata,MPI_SUM)); /* rank 1->0, add leafdata to rootdata */
51*9566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd(sf,MPI_SIGNED_CHAR,leafdata,rootdata,MPI_SUM));
5281b34b04SJunchao Zhang   if (!rank) {
532c71b3e2SJacob Faibussowitsch     PetscCheckFalse(rootdata[0] != 22 || rootdata[nroots-1] != 24,PETSC_COMM_SELF,PETSC_ERR_PLIB,"SF: wrong results");
5481b34b04SJunchao Zhang   }
5581b34b04SJunchao Zhang 
56*9566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rootdata,leafdata));
57*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(iremote));
58*9566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sf));
5981b34b04SJunchao Zhang 
6081b34b04SJunchao Zhang   /* Test VecScatter */
61*9566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
62*9566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&y));
63*9566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x,rank==0? n : 64,PETSC_DECIDE));
64*9566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(y,rank==0? 64 : n,PETSC_DECIDE));
65*9566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
66*9566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(y));
6781b34b04SJunchao Zhang 
68*9566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(x,&rstart,&rend));
69*9566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF,rend-rstart,rstart,1,&ix));
70*9566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x,ix,y,ix,&vscat));
7181b34b04SJunchao Zhang 
72*9566063dSJacob Faibussowitsch   PetscCall(VecSet(x,3.0));
73*9566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
74*9566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
7581b34b04SJunchao Zhang 
76*9566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(vscat,y,x,ADD_VALUES,SCATTER_REVERSE));
77*9566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(vscat,y,x,ADD_VALUES,SCATTER_REVERSE));
7881b34b04SJunchao Zhang 
79*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x,&xv));
802c71b3e2SJacob Faibussowitsch   PetscCheckFalse(xv[0] != 6.0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"VecScatter: wrong results");
81*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x,&xv));
82*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
83*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
84*9566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&vscat));
85*9566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ix));
8681b34b04SJunchao Zhang 
87*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
88b122ec5aSJacob Faibussowitsch   return 0;
8981b34b04SJunchao Zhang }
9081b34b04SJunchao Zhang 
9181b34b04SJunchao Zhang /**TEST
9281b34b04SJunchao Zhang    test:
9381b34b04SJunchao Zhang      requires: defined(PETSC_HAVE_MPI_LARGE_COUNT) defined(PETSC_USE_64BIT_INDICES)
9481b34b04SJunchao Zhang      TODO: need a machine with big memory (~150GB) to run the test
9581b34b04SJunchao Zhang      nsize: 2
9681b34b04SJunchao Zhang      args: -sf_type {{basic neighbor}}
9781b34b04SJunchao Zhang 
9881b34b04SJunchao Zhang TEST**/
99