xref: /petsc/src/vec/is/sf/tests/ex17.c (revision c5853193be7e0cd23913bbf6af39be2e963e2513)
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 
209566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
2308401ef6SPierre Jolivet   PetscCheck(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 */
269566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PETSC_COMM_WORLD,&sf));
279566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sf));
2881b34b04SJunchao Zhang 
29*c5853193SPierre Jolivet   if (rank == 0) {
3081b34b04SJunchao Zhang     nroots  = n;
3181b34b04SJunchao Zhang     nleaves = 0;
3281b34b04SJunchao Zhang   } else {
3381b34b04SJunchao Zhang     nroots  = 0;
3481b34b04SJunchao Zhang     nleaves = n;
359566063dSJacob 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   }
419566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(sf,nroots,nleaves,NULL,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES));
429566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nroots,&rootdata,nleaves,&leafdata));
43*c5853193SPierre Jolivet   if (rank == 0) {
4481b34b04SJunchao Zhang     memset(rootdata,11,nroots);
4581b34b04SJunchao Zhang     rootdata[nroots-1] = 12; /* Use a different value at the end */
4681b34b04SJunchao Zhang   }
4781b34b04SJunchao Zhang 
489566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sf,MPI_SIGNED_CHAR,rootdata,leafdata,MPI_REPLACE)); /* rank 0->1, bcast rootdata to leafdata */
499566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sf,MPI_SIGNED_CHAR,rootdata,leafdata,MPI_REPLACE));
509566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(sf,MPI_SIGNED_CHAR,leafdata,rootdata,MPI_SUM)); /* rank 1->0, add leafdata to rootdata */
519566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd(sf,MPI_SIGNED_CHAR,leafdata,rootdata,MPI_SUM));
52*c5853193SPierre Jolivet   PetscCheck(rank != 0 || (rootdata[0] == 22 && rootdata[nroots-1] == 24),PETSC_COMM_SELF,PETSC_ERR_PLIB,"SF: wrong results");
5381b34b04SJunchao Zhang 
549566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rootdata,leafdata));
559566063dSJacob Faibussowitsch   PetscCall(PetscFree(iremote));
569566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sf));
5781b34b04SJunchao Zhang 
5881b34b04SJunchao Zhang   /* Test VecScatter */
599566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
609566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&y));
619566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x,rank==0? n : 64,PETSC_DECIDE));
629566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(y,rank==0? 64 : n,PETSC_DECIDE));
639566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
649566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(y));
6581b34b04SJunchao Zhang 
669566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(x,&rstart,&rend));
679566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF,rend-rstart,rstart,1,&ix));
689566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x,ix,y,ix,&vscat));
6981b34b04SJunchao Zhang 
709566063dSJacob Faibussowitsch   PetscCall(VecSet(x,3.0));
719566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
729566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
7381b34b04SJunchao Zhang 
749566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(vscat,y,x,ADD_VALUES,SCATTER_REVERSE));
759566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(vscat,y,x,ADD_VALUES,SCATTER_REVERSE));
7681b34b04SJunchao Zhang 
779566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x,&xv));
7808401ef6SPierre Jolivet   PetscCheck(xv[0] == 6.0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"VecScatter: wrong results");
799566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x,&xv));
809566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
829566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&vscat));
839566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ix));
8481b34b04SJunchao Zhang 
859566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
86b122ec5aSJacob Faibussowitsch   return 0;
8781b34b04SJunchao Zhang }
8881b34b04SJunchao Zhang 
8981b34b04SJunchao Zhang /**TEST
9081b34b04SJunchao Zhang    test:
9181b34b04SJunchao Zhang      requires: defined(PETSC_HAVE_MPI_LARGE_COUNT) defined(PETSC_USE_64BIT_INDICES)
9281b34b04SJunchao Zhang      TODO: need a machine with big memory (~150GB) to run the test
9381b34b04SJunchao Zhang      nsize: 2
9481b34b04SJunchao Zhang      args: -sf_type {{basic neighbor}}
9581b34b04SJunchao Zhang 
9681b34b04SJunchao Zhang TEST**/
97