xref: /petsc/src/vec/is/sf/tests/ex17.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode    ierr;
981b34b04SJunchao Zhang   PetscSF           sf;
1081b34b04SJunchao Zhang   PetscInt          i,nroots,nleaves;
1181b34b04SJunchao Zhang   PetscInt          n = (1ULL<<31) + 1024; /* a little over 2G elements */
1281b34b04SJunchao Zhang   PetscSFNode       *iremote = NULL;
1381b34b04SJunchao Zhang   PetscMPIInt       rank,size;
1481b34b04SJunchao Zhang   char              *rootdata=NULL,*leafdata=NULL;
1581b34b04SJunchao Zhang   Vec               x,y;
1681b34b04SJunchao Zhang   VecScatter        vscat;
1781b34b04SJunchao Zhang   PetscInt          rstart,rend;
1881b34b04SJunchao Zhang   IS                ix;
1981b34b04SJunchao Zhang   const PetscScalar *xv;
2081b34b04SJunchao Zhang 
2181b34b04SJunchao Zhang   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
22*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
23*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
242c71b3e2SJacob Faibussowitsch   PetscCheckFalse(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 */
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFCreate(PETSC_COMM_WORLD,&sf));
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetFromOptions(sf));
2981b34b04SJunchao Zhang 
3081b34b04SJunchao Zhang   if (!rank) {
3181b34b04SJunchao Zhang     nroots  = n;
3281b34b04SJunchao Zhang     nleaves = 0;
3381b34b04SJunchao Zhang   } else {
3481b34b04SJunchao Zhang     nroots  = 0;
3581b34b04SJunchao Zhang     nleaves = n;
36*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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   }
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetGraph(sf,nroots,nleaves,NULL,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES));
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(nroots,&rootdata,nleaves,&leafdata));
4481b34b04SJunchao Zhang   if (!rank) {
4581b34b04SJunchao Zhang     memset(rootdata,11,nroots);
4681b34b04SJunchao Zhang     rootdata[nroots-1] = 12; /* Use a different value at the end */
4781b34b04SJunchao Zhang   }
4881b34b04SJunchao Zhang 
49*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastBegin(sf,MPI_SIGNED_CHAR,rootdata,leafdata,MPI_REPLACE)); /* rank 0->1, bcast rootdata to leafdata */
50*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastEnd(sf,MPI_SIGNED_CHAR,rootdata,leafdata,MPI_REPLACE));
51*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFReduceBegin(sf,MPI_SIGNED_CHAR,leafdata,rootdata,MPI_SUM)); /* rank 1->0, add leafdata to rootdata */
52*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFReduceEnd(sf,MPI_SIGNED_CHAR,leafdata,rootdata,MPI_SUM));
5381b34b04SJunchao Zhang   if (!rank) {
542c71b3e2SJacob Faibussowitsch     PetscCheckFalse(rootdata[0] != 22 || rootdata[nroots-1] != 24,PETSC_COMM_SELF,PETSC_ERR_PLIB,"SF: wrong results");
5581b34b04SJunchao Zhang   }
5681b34b04SJunchao Zhang 
57*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(rootdata,leafdata));
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(iremote));
59*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFDestroy(&sf));
6081b34b04SJunchao Zhang 
6181b34b04SJunchao Zhang   /* Test VecScatter */
62*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x));
63*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&y));
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(x,rank==0? n : 64,PETSC_DECIDE));
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(y,rank==0? 64 : n,PETSC_DECIDE));
66*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(x));
67*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(y));
6881b34b04SJunchao Zhang 
69*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetOwnershipRange(x,&rstart,&rend));
70*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateStride(PETSC_COMM_SELF,rend-rstart,rstart,1,&ix));
71*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(x,ix,y,ix,&vscat));
7281b34b04SJunchao Zhang 
73*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(x,3.0));
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
75*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,x,y,INSERT_VALUES,SCATTER_FORWARD));
7681b34b04SJunchao Zhang 
77*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(vscat,y,x,ADD_VALUES,SCATTER_REVERSE));
78*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(vscat,y,x,ADD_VALUES,SCATTER_REVERSE));
7981b34b04SJunchao Zhang 
80*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(x,&xv));
812c71b3e2SJacob Faibussowitsch   PetscCheckFalse(xv[0] != 6.0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"VecScatter: wrong results");
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(x,&xv));
83*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
84*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
85*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&vscat));
86*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&ix));
8781b34b04SJunchao Zhang 
8881b34b04SJunchao Zhang   ierr = PetscFinalize();
8981b34b04SJunchao Zhang   return ierr;
9081b34b04SJunchao Zhang }
9181b34b04SJunchao Zhang 
9281b34b04SJunchao Zhang /**TEST
9381b34b04SJunchao Zhang    test:
9481b34b04SJunchao Zhang      requires: defined(PETSC_HAVE_MPI_LARGE_COUNT) defined(PETSC_USE_64BIT_INDICES)
9581b34b04SJunchao Zhang      TODO: need a machine with big memory (~150GB) to run the test
9681b34b04SJunchao Zhang      nsize: 2
9781b34b04SJunchao Zhang      args: -sf_type {{basic neighbor}}
9881b34b04SJunchao Zhang 
9981b34b04SJunchao Zhang TEST**/
100