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