1*c4762a1bSJed Brown static const char help[] = "Test overlapped communication on a single star forest (PetscSF)\n\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscvec.h> 4*c4762a1bSJed Brown #include <petscsf.h> 5*c4762a1bSJed Brown #include <petscviewer.h> 6*c4762a1bSJed Brown 7*c4762a1bSJed Brown int main(int argc, char **argv) 8*c4762a1bSJed Brown { 9*c4762a1bSJed Brown PetscInt ierr; 10*c4762a1bSJed Brown PetscSF sf; 11*c4762a1bSJed Brown Vec A,Aout; 12*c4762a1bSJed Brown Vec B,Bout; 13*c4762a1bSJed Brown PetscScalar *bufA; 14*c4762a1bSJed Brown PetscScalar *bufAout; 15*c4762a1bSJed Brown PetscScalar *bufB; 16*c4762a1bSJed Brown PetscScalar *bufBout; 17*c4762a1bSJed Brown PetscMPIInt rank, size; 18*c4762a1bSJed Brown PetscInt nroots, nleaves; 19*c4762a1bSJed Brown PetscInt i; 20*c4762a1bSJed Brown PetscInt *ilocal; 21*c4762a1bSJed Brown PetscSFNode *iremote; 22*c4762a1bSJed Brown 23*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 24*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 25*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 26*c4762a1bSJed Brown 27*c4762a1bSJed Brown if (size != 2) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_USER, "Only coded for two MPI processes\n"); 28*c4762a1bSJed Brown 29*c4762a1bSJed Brown ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); 30*c4762a1bSJed Brown ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown nleaves = 2; 33*c4762a1bSJed Brown nroots = 1; 34*c4762a1bSJed Brown ierr = PetscMalloc1(nleaves,&ilocal);CHKERRQ(ierr); 35*c4762a1bSJed Brown 36*c4762a1bSJed Brown for (i = 0; i<nleaves; i++) { 37*c4762a1bSJed Brown ilocal[i] = i; 38*c4762a1bSJed Brown } 39*c4762a1bSJed Brown 40*c4762a1bSJed Brown ierr = PetscMalloc1(nleaves,&iremote);CHKERRQ(ierr); 41*c4762a1bSJed Brown if (rank == 0) { 42*c4762a1bSJed Brown iremote[0].rank = 0; 43*c4762a1bSJed Brown iremote[0].index = 0; 44*c4762a1bSJed Brown iremote[1].rank = 1; 45*c4762a1bSJed Brown iremote[1].index = 0; 46*c4762a1bSJed Brown } else { 47*c4762a1bSJed Brown iremote[0].rank = 1; 48*c4762a1bSJed Brown iremote[0].index = 0; 49*c4762a1bSJed Brown iremote[1].rank = 0; 50*c4762a1bSJed Brown iremote[1].index = 0; 51*c4762a1bSJed Brown } 52*c4762a1bSJed Brown ierr = PetscSFSetGraph(sf,nroots,nleaves,ilocal,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 53*c4762a1bSJed Brown ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 54*c4762a1bSJed Brown ierr = PetscSFView(sf,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 55*c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 56*c4762a1bSJed Brown ierr = VecSetSizes(A,2,PETSC_DETERMINE);CHKERRQ(ierr); 57*c4762a1bSJed Brown ierr = VecSetFromOptions(A);CHKERRQ(ierr); 58*c4762a1bSJed Brown ierr = VecSetUp(A);CHKERRQ(ierr); 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown ierr = VecDuplicate(A,&B);CHKERRQ(ierr); 61*c4762a1bSJed Brown ierr = VecDuplicate(A,&Aout);CHKERRQ(ierr); 62*c4762a1bSJed Brown ierr = VecDuplicate(A,&Bout);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = VecGetArray(A,&bufA);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = VecGetArray(B,&bufB);CHKERRQ(ierr); 65*c4762a1bSJed Brown for (i=0; i<2; i++) { 66*c4762a1bSJed Brown bufA[i] = (PetscScalar)rank; 67*c4762a1bSJed Brown bufB[i] = (PetscScalar)(rank) + 10.0; 68*c4762a1bSJed Brown } 69*c4762a1bSJed Brown ierr = VecRestoreArray(A,&bufA);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = VecRestoreArray(B,&bufB);CHKERRQ(ierr); 71*c4762a1bSJed Brown 72*c4762a1bSJed Brown ierr = VecGetArrayRead(A,(const PetscScalar**)&bufA);CHKERRQ(ierr); 73*c4762a1bSJed Brown ierr = VecGetArrayRead(B,(const PetscScalar**)&bufB);CHKERRQ(ierr); 74*c4762a1bSJed Brown ierr = VecGetArray(Aout,&bufAout);CHKERRQ(ierr); 75*c4762a1bSJed Brown ierr = VecGetArray(Bout,&bufBout);CHKERRQ(ierr); 76*c4762a1bSJed Brown ierr = PetscSFBcastBegin(sf,MPIU_SCALAR,(const void*)bufA,(void *)bufAout);CHKERRQ(ierr); 77*c4762a1bSJed Brown ierr = PetscSFBcastBegin(sf,MPIU_SCALAR,(const void*)bufB,(void *)bufBout);CHKERRQ(ierr); 78*c4762a1bSJed Brown ierr = PetscSFBcastEnd(sf,MPIU_SCALAR,(const void*)bufA,(void *)bufAout);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = PetscSFBcastEnd(sf,MPIU_SCALAR,(const void*)bufB,(void *)bufBout);CHKERRQ(ierr); 80*c4762a1bSJed Brown ierr = VecRestoreArrayRead(A,(const PetscScalar**)&bufA);CHKERRQ(ierr); 81*c4762a1bSJed Brown ierr = VecRestoreArrayRead(B,(const PetscScalar**)&bufB);CHKERRQ(ierr); 82*c4762a1bSJed Brown ierr = VecRestoreArray(Aout,&bufAout);CHKERRQ(ierr); 83*c4762a1bSJed Brown ierr = VecRestoreArray(Bout,&bufBout);CHKERRQ(ierr); 84*c4762a1bSJed Brown 85*c4762a1bSJed Brown ierr = VecView(Aout,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 86*c4762a1bSJed Brown ierr = VecView(Bout,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 87*c4762a1bSJed Brown ierr = VecDestroy(&A);CHKERRQ(ierr); 88*c4762a1bSJed Brown ierr = VecDestroy(&B);CHKERRQ(ierr); 89*c4762a1bSJed Brown ierr = VecDestroy(&Aout);CHKERRQ(ierr); 90*c4762a1bSJed Brown ierr = VecDestroy(&Bout);CHKERRQ(ierr); 91*c4762a1bSJed Brown ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 92*c4762a1bSJed Brown 93*c4762a1bSJed Brown ierr = PetscFinalize(); 94*c4762a1bSJed Brown return ierr; 95*c4762a1bSJed Brown } 96*c4762a1bSJed Brown 97*c4762a1bSJed Brown 98*c4762a1bSJed Brown /*TEST 99*c4762a1bSJed Brown 100*c4762a1bSJed Brown test: 101*c4762a1bSJed Brown suffix: basic 102*c4762a1bSJed Brown nsize: 2 103*c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 104*c4762a1bSJed Brown args: -sf_type basic 105*c4762a1bSJed Brown 106*c4762a1bSJed Brown test: 107*c4762a1bSJed Brown suffix: window 108*c4762a1bSJed Brown nsize: 2 109*c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 110*c4762a1bSJed Brown output_file: output/ex2_basic.out 111*c4762a1bSJed Brown args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 112*c4762a1bSJed Brown requires: define(PETSC_HAVE_MPI_ONE_SIDED) 113*c4762a1bSJed Brown 114*c4762a1bSJed Brown # The nightly test suite with MPICH uses ch3:sock, which is broken when winsize == 0 in some of the processes 115*c4762a1bSJed Brown test: 116*c4762a1bSJed Brown suffix: window_shared 117*c4762a1bSJed Brown nsize: 2 118*c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 119*c4762a1bSJed Brown output_file: output/ex2_basic.out 120*c4762a1bSJed Brown args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor shared 121*c4762a1bSJed Brown requires: define(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) !define(PETSC_HAVE_MPICH_NUMVERSION) define(PETSC_HAVE_MPI_ONE_SIDED) 122*c4762a1bSJed Brown 123*c4762a1bSJed Brown TEST*/ 124