Actual source code: ex2.c

  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>

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

 22:   PetscInitialize(&argc,&argv,NULL,help);
 23:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 24:   MPI_Comm_size(PETSC_COMM_WORLD,&size);


 28:   PetscSFCreate(PETSC_COMM_WORLD,&sf);
 29:   PetscSFSetFromOptions(sf);

 31:   nleaves = 2;
 32:   nroots = 1;
 33:   PetscMalloc1(nleaves,&ilocal);

 35:   for (i = 0; i<nleaves; i++) {
 36:     ilocal[i] = i;
 37:   }

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

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

 71:   VecGetArrayRead(A,(const PetscScalar**)&bufA);
 72:   VecGetArrayRead(B,(const PetscScalar**)&bufB);
 73:   VecGetArray(Aout,&bufAout);
 74:   VecGetArray(Bout,&bufBout);
 75:   PetscSFBcastBegin(sf,MPIU_SCALAR,(const void*)bufA,(void *)bufAout,MPI_REPLACE);
 76:   PetscSFBcastBegin(sf,MPIU_SCALAR,(const void*)bufB,(void *)bufBout,MPI_REPLACE);
 77:   PetscSFBcastEnd(sf,MPIU_SCALAR,(const void*)bufA,(void *)bufAout,MPI_REPLACE);
 78:   PetscSFBcastEnd(sf,MPIU_SCALAR,(const void*)bufB,(void *)bufBout,MPI_REPLACE);
 79:   VecRestoreArrayRead(A,(const PetscScalar**)&bufA);
 80:   VecRestoreArrayRead(B,(const PetscScalar**)&bufB);
 81:   VecRestoreArray(Aout,&bufAout);
 82:   VecRestoreArray(Bout,&bufBout);

 84:   VecView(Aout,PETSC_VIEWER_STDOUT_WORLD);
 85:   VecView(Bout,PETSC_VIEWER_STDOUT_WORLD);
 86:   VecDestroy(&A);
 87:   VecDestroy(&B);
 88:   VecDestroy(&Aout);
 89:   VecDestroy(&Bout);
 90:   PetscSFDestroy(&sf);

 92:   PetscFinalize();
 93:   return 0;
 94: }

 96: /*TEST

 98:    test:
 99:       suffix: basic
100:       nsize: 2
101:       filter: grep -v "type" | grep -v "sort"
102:       args: -sf_type basic

104:    test:
105:       suffix: window
106:       nsize: 2
107:       filter: grep -v "type" | grep -v "sort"
108:       output_file: output/ex2_basic.out
109:       args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
110:       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)

112:    # The nightly test suite with MPICH uses ch3:sock, which is broken when winsize == 0 in some of the processes
113:    test:
114:       suffix: window_shared
115:       nsize: 2
116:       filter: grep -v "type" | grep -v "sort"
117:       output_file: output/ex2_basic.out
118:       args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor shared
119:       requires: defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) !defined(PETSC_HAVE_MPICH_NUMVERSION) defined(PETSC_HAVE_MPI_ONE_SIDED)

121: TEST*/