Actual source code: ex2.c

petsc-3.11.4 2019-09-28
Report Typos and Errors
  1: static const char help[] = "Test embedded sf with one leaf data item connnected to multiple roots\n\n";

  3:  #include <petscsf.h>

  5: int main(int argc,char **argv)
  6: {
  7:   PetscSF        sf,newsf;
  8:   PetscInt       i,nroots,nleaves,ilocal[2],leafdata,rootdata[2],nselected,selected,errors=0;
  9:   PetscSFNode    iremote[2];
 10:   PetscMPIInt    myrank,next,nproc;

 13:   PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
 14:   MPI_Comm_size(PETSC_COMM_WORLD,&nproc);
 15:   MPI_Comm_rank(PETSC_COMM_WORLD,&myrank);

 17:   /* Create an SF that each process has two roots and two leaves. The two leaves are connected to the two
 18:      roots on next process. The special thing is each process only has one data item in its leaf data buffer.
 19:    */
 20:   nroots    = 2;
 21:   nleaves   = 2;
 22:   ilocal[0] = 0; /* One leaf data item serves two leaves */
 23:   ilocal[1] = 0;
 24:   next      = (myrank+1)%nproc;
 25:   for (i=0; i<nleaves; i++) {
 26:     iremote[i].rank  = next;
 27:     iremote[i].index = i;
 28:   }

 30:   PetscSFCreate(PETSC_COMM_WORLD,&sf);
 31:   PetscSFSetFromOptions(sf);
 32:   PetscSFSetGraph(sf,nroots,nleaves,ilocal,PETSC_COPY_VALUES,iremote,PETSC_COPY_VALUES);

 34:   /* Create an embedded sf by only selecting the first root on each process */
 35:   nselected = 1;
 36:   selected  = 0;
 37:   PetscSFCreateEmbeddedSF(sf,nselected,&selected,&newsf);

 39:   /* Do reduce */
 40:   leafdata    = 1;
 41:   rootdata[0] = rootdata[1] = 0;
 42:   PetscSFReduceBegin(newsf,MPIU_INT,&leafdata,rootdata,MPIU_REPLACE);
 43:   PetscSFReduceEnd(newsf,MPIU_INT,&leafdata,rootdata,MPIU_REPLACE);

 45:   /* Check rootdata */
 46:   if (rootdata[0] != 1 || rootdata[1] != 0) errors = 1;
 47:   MPI_Allreduce(MPI_IN_PLACE,&errors,1,MPIU_INT,MPI_SUM,PETSC_COMM_WORLD);
 48:   if (errors) {PetscPrintf(PETSC_COMM_WORLD,"Error: Unexpected rootdata on processors\n");}

 50:   PetscSFDestroy(&sf);
 51:   PetscSFDestroy(&newsf);
 52:   PetscFinalize();
 53:   return ierr;
 54: }

 56: /*TEST

 58:    test:
 59:       nsize: 2
 60: TEST*/