Actual source code: ex2.c

petsc-3.7.7 2017-09-25
Report Typos and Errors
  1: static const char help[] = "Test overlapped communication on a single star forest (PetscSF)\n\n";

  3: #include <petscvec.h>
  4: #include <petscsf.h>
  5: #include <petscviewer.h>

  9: int main(int argc, char **argv)
 10: {
 11:   PetscInt    ierr;
 12:   PetscSF     sf;
 13:   Vec         A,Aout;
 14:   Vec         B,Bout;
 15:   PetscScalar *bufA;
 16:   PetscScalar *bufAout;
 17:   PetscScalar *bufB;
 18:   PetscScalar *bufBout;
 19:   PetscMPIInt rank, size;
 20:   PetscInt    nroots, nleaves;
 21:   PetscInt    i;
 22:   PetscInt    *ilocal;
 23:   PetscSFNode *iremote;

 25:   PetscInitialize(&argc,&argv,NULL,help);

 27:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 28:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 30:   if (size != 2) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_USER, "Only coded for two MPI processes\n");

 32:   PetscSFCreate(PETSC_COMM_WORLD,&sf);
 33:   PetscSFSetFromOptions(sf);

 35:   nleaves = 2;
 36:   nroots = 1;
 37:   PetscMalloc1(nleaves,&ilocal);

 39:   for (i = 0; i<nleaves; i++) {
 40:     ilocal[i] = i;
 41:   }

 43:   PetscMalloc1(nleaves,&iremote);
 44:   if (rank == 0) {
 45:     iremote[0].rank = 0;
 46:     iremote[0].index = 0;
 47:     iremote[1].rank = 1;
 48:     iremote[1].index = 0;
 49:   } else {
 50:     iremote[0].rank = 1;
 51:     iremote[0].index = 0;
 52:     iremote[1].rank = 0;
 53:     iremote[1].index = 0;
 54:   }
 55:   PetscSFSetGraph(sf,nroots,nleaves,ilocal,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);
 56:   PetscSFSetUp(sf);
 57:   PetscSFView(sf,PETSC_VIEWER_STDOUT_WORLD);
 58:   VecCreate(PETSC_COMM_WORLD,&A);
 59:   VecSetSizes(A,2,PETSC_DETERMINE);
 60:   VecSetFromOptions(A);
 61:   VecSetUp(A);

 63:   VecDuplicate(A,&B);
 64:   VecDuplicate(A,&Aout);
 65:   VecDuplicate(A,&Bout);
 66:   VecGetArray(A,&bufA);
 67:   VecGetArray(B,&bufB);
 68:   for (i=0; i<2; i++) {
 69:     bufA[i] = (PetscScalar)rank;
 70:     bufB[i] = (PetscScalar)(rank) + 10.0;
 71:   }
 72:   VecRestoreArray(A,&bufA);
 73:   VecRestoreArray(B,&bufB);

 75:   VecGetArrayRead(A,(const PetscScalar**)&bufA);
 76:   VecGetArrayRead(B,(const PetscScalar**)&bufB);
 77:   VecGetArray(Aout,&bufAout);
 78:   VecGetArray(Bout,&bufBout);
 79:   PetscSFBcastBegin(sf,MPIU_SCALAR,(const void*)bufA,(void *)bufAout);
 80:   PetscSFBcastBegin(sf,MPIU_SCALAR,(const void*)bufB,(void *)bufBout);
 81:   PetscSFBcastEnd(sf,MPIU_SCALAR,(const void*)bufA,(void *)bufAout);
 82:   PetscSFBcastEnd(sf,MPIU_SCALAR,(const void*)bufB,(void *)bufBout);
 83:   VecRestoreArrayRead(A,(const PetscScalar**)&bufA);
 84:   VecRestoreArrayRead(B,(const PetscScalar**)&bufB);
 85:   VecRestoreArray(Aout,&bufAout);
 86:   VecRestoreArray(Bout,&bufBout);

 88:   VecView(Aout,PETSC_VIEWER_STDOUT_WORLD);
 89:   VecView(Bout,PETSC_VIEWER_STDOUT_WORLD);
 90:   VecDestroy(&A);
 91:   VecDestroy(&B);
 92:   VecDestroy(&Aout);
 93:   VecDestroy(&Bout);
 94:   PetscSFDestroy(&sf);

 96:   PetscFinalize();
 97:   return 0;
 98: }