1*33e3ecb2SToby Isaac const char help[] = "Test overlapping PetscSF communication with empty roots and leaves"; 2*33e3ecb2SToby Isaac 3*33e3ecb2SToby Isaac #include <petscsf.h> 4*33e3ecb2SToby Isaac 5*33e3ecb2SToby Isaac static PetscErrorCode testOverlappingCommunication(PetscSF sf) 6*33e3ecb2SToby Isaac { 7*33e3ecb2SToby Isaac PetscInt nroots, maxleaf; 8*33e3ecb2SToby Isaac PetscInt *leafa, *leafb, *roota, *rootb; 9*33e3ecb2SToby Isaac 10*33e3ecb2SToby Isaac PetscFunctionBegin; 11*33e3ecb2SToby Isaac PetscCall(PetscSFSetUp(sf)); 12*33e3ecb2SToby Isaac PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL)); 13*33e3ecb2SToby Isaac PetscCall(PetscSFGetLeafRange(sf, NULL, &maxleaf)); 14*33e3ecb2SToby Isaac PetscCall(PetscMalloc4(nroots, &roota, nroots, &rootb, maxleaf + 1, &leafa, maxleaf + 1, &leafb)); 15*33e3ecb2SToby Isaac 16*33e3ecb2SToby Isaac // test reduce 17*33e3ecb2SToby Isaac for (PetscInt i = 0; i < nroots; i++) roota[i] = 0; 18*33e3ecb2SToby Isaac for (PetscInt i = 0; i < nroots; i++) rootb[i] = 0; 19*33e3ecb2SToby Isaac for (PetscInt i = 0; i < maxleaf + 1; i++) leafa[i] = (i + 1); 20*33e3ecb2SToby Isaac for (PetscInt i = 0; i < maxleaf + 1; i++) leafb[i] = -(i + 1); 21*33e3ecb2SToby Isaac 22*33e3ecb2SToby Isaac PetscCall(PetscSFReduceBegin(sf, MPIU_INT, leafa, roota, MPI_REPLACE)); 23*33e3ecb2SToby Isaac PetscCall(PetscSFReduceBegin(sf, MPIU_INT, leafb, rootb, MPI_REPLACE)); 24*33e3ecb2SToby Isaac PetscCall(PetscSFReduceEnd(sf, MPIU_INT, leafa, roota, MPI_REPLACE)); 25*33e3ecb2SToby Isaac PetscCall(PetscSFReduceEnd(sf, MPIU_INT, leafb, rootb, MPI_REPLACE)); 26*33e3ecb2SToby Isaac for (PetscInt i = 0; i < nroots; i++) PetscCheck(roota[i] >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscSFReduce in (A,B,A,B) order crosses separate reductions"); 27*33e3ecb2SToby Isaac for (PetscInt i = 0; i < nroots; i++) PetscCheck(rootb[i] <= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscSFReduce in (A,B,A,B) order crosses separate reductions"); 28*33e3ecb2SToby Isaac 29*33e3ecb2SToby Isaac for (PetscInt i = 0; i < nroots; i++) roota[i] = 0; 30*33e3ecb2SToby Isaac for (PetscInt i = 0; i < nroots; i++) rootb[i] = 0; 31*33e3ecb2SToby Isaac 32*33e3ecb2SToby Isaac PetscCall(PetscSFReduceBegin(sf, MPIU_INT, leafa, roota, MPI_REPLACE)); 33*33e3ecb2SToby Isaac PetscCall(PetscSFReduceBegin(sf, MPIU_INT, leafb, rootb, MPI_REPLACE)); 34*33e3ecb2SToby Isaac PetscCall(PetscSFReduceEnd(sf, MPIU_INT, leafb, rootb, MPI_REPLACE)); 35*33e3ecb2SToby Isaac PetscCall(PetscSFReduceEnd(sf, MPIU_INT, leafa, roota, MPI_REPLACE)); 36*33e3ecb2SToby Isaac 37*33e3ecb2SToby Isaac for (PetscInt i = 0; i < nroots; i++) PetscCheck(roota[i] >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscSFReduce in (A,B,B,A) order crosses separate reductions"); 38*33e3ecb2SToby Isaac for (PetscInt i = 0; i < nroots; i++) PetscCheck(rootb[i] <= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscSFReduce in (A,B,B,A) order crosses separate reductions"); 39*33e3ecb2SToby Isaac 40*33e3ecb2SToby Isaac // test bcast 41*33e3ecb2SToby Isaac for (PetscInt i = 0; i < nroots; i++) roota[i] = (i + 1); 42*33e3ecb2SToby Isaac for (PetscInt i = 0; i < nroots; i++) rootb[i] = -(i + 1); 43*33e3ecb2SToby Isaac for (PetscInt i = 0; i < maxleaf + 1; i++) leafa[i] = 0; 44*33e3ecb2SToby Isaac for (PetscInt i = 0; i < maxleaf + 1; i++) leafb[i] = 0; 45*33e3ecb2SToby Isaac 46*33e3ecb2SToby Isaac PetscCall(PetscSFBcastBegin(sf, MPIU_INT, roota, leafa, MPI_REPLACE)); 47*33e3ecb2SToby Isaac PetscCall(PetscSFBcastBegin(sf, MPIU_INT, rootb, leafb, MPI_REPLACE)); 48*33e3ecb2SToby Isaac PetscCall(PetscSFBcastEnd(sf, MPIU_INT, roota, leafa, MPI_REPLACE)); 49*33e3ecb2SToby Isaac PetscCall(PetscSFBcastEnd(sf, MPIU_INT, rootb, leafb, MPI_REPLACE)); 50*33e3ecb2SToby Isaac 51*33e3ecb2SToby Isaac for (PetscInt i = 0; i < maxleaf + 1; i++) PetscCheck(leafa[i] >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscSFBcast in (A,B,A,B) order crosses separate broadcasts"); 52*33e3ecb2SToby Isaac for (PetscInt i = 0; i < maxleaf + 1; i++) PetscCheck(leafb[i] <= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscSFBcast in (A,B,A,B) order crosses separate broadcasts"); 53*33e3ecb2SToby Isaac 54*33e3ecb2SToby Isaac for (PetscInt i = 0; i < maxleaf + 1; i++) leafa[i] = 0; 55*33e3ecb2SToby Isaac for (PetscInt i = 0; i < maxleaf + 1; i++) leafb[i] = 0; 56*33e3ecb2SToby Isaac 57*33e3ecb2SToby Isaac PetscCall(PetscSFBcastBegin(sf, MPIU_INT, roota, leafa, MPI_REPLACE)); 58*33e3ecb2SToby Isaac PetscCall(PetscSFBcastBegin(sf, MPIU_INT, rootb, leafb, MPI_REPLACE)); 59*33e3ecb2SToby Isaac PetscCall(PetscSFBcastEnd(sf, MPIU_INT, rootb, leafb, MPI_REPLACE)); 60*33e3ecb2SToby Isaac PetscCall(PetscSFBcastEnd(sf, MPIU_INT, roota, leafa, MPI_REPLACE)); 61*33e3ecb2SToby Isaac 62*33e3ecb2SToby Isaac for (PetscInt i = 0; i < maxleaf + 1; i++) PetscCheck(leafa[i] >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscSFBcast in (A,B,B,A) order crosses separate broadcasts"); 63*33e3ecb2SToby Isaac for (PetscInt i = 0; i < maxleaf + 1; i++) PetscCheck(leafb[i] <= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscSFBcast in (A,B,B,A) order crosses separate broadcasts"); 64*33e3ecb2SToby Isaac 65*33e3ecb2SToby Isaac PetscCall(PetscFree4(roota, rootb, leafa, leafb)); 66*33e3ecb2SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 67*33e3ecb2SToby Isaac } 68*33e3ecb2SToby Isaac 69*33e3ecb2SToby Isaac static PetscErrorCode createSparseSF(MPI_Comm comm, PetscSF *sf) 70*33e3ecb2SToby Isaac { 71*33e3ecb2SToby Isaac PetscMPIInt rank; 72*33e3ecb2SToby Isaac PetscInt nroots, nleaves; 73*33e3ecb2SToby Isaac PetscSFNode remote; 74*33e3ecb2SToby Isaac 75*33e3ecb2SToby Isaac PetscFunctionBegin; 76*33e3ecb2SToby Isaac PetscCallMPI(MPI_Comm_rank(comm, &rank)); 77*33e3ecb2SToby Isaac PetscCall(PetscSFCreate(comm, sf)); 78*33e3ecb2SToby Isaac nroots = (rank & 1) ? 1 : 0; 79*33e3ecb2SToby Isaac nleaves = (rank & 2) ? 1 : 0; 80*33e3ecb2SToby Isaac remote.rank = -1; 81*33e3ecb2SToby Isaac remote.index = 0; 82*33e3ecb2SToby Isaac if (nleaves == 1) { remote.rank = (rank & 1) ? (rank ^ 2) : (rank ^ 1); } 83*33e3ecb2SToby Isaac PetscCall(PetscSFSetGraph(*sf, nroots, nleaves, NULL, PETSC_COPY_VALUES, &remote, PETSC_COPY_VALUES)); 84*33e3ecb2SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 85*33e3ecb2SToby Isaac } 86*33e3ecb2SToby Isaac 87*33e3ecb2SToby Isaac int main(int argc, char **argv) 88*33e3ecb2SToby Isaac { 89*33e3ecb2SToby Isaac PetscSF sf; 90*33e3ecb2SToby Isaac 91*33e3ecb2SToby Isaac PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 92*33e3ecb2SToby Isaac PetscCall(createSparseSF(PETSC_COMM_WORLD, &sf)); 93*33e3ecb2SToby Isaac PetscCall(PetscSFSetFromOptions(sf)); 94*33e3ecb2SToby Isaac PetscCall(testOverlappingCommunication(sf)); 95*33e3ecb2SToby Isaac PetscCall(PetscSFDestroy(&sf)); 96*33e3ecb2SToby Isaac PetscCall(PetscFinalize()); 97*33e3ecb2SToby Isaac return 0; 98*33e3ecb2SToby Isaac } 99*33e3ecb2SToby Isaac 100*33e3ecb2SToby Isaac /*TEST 101*33e3ecb2SToby Isaac 102*33e3ecb2SToby Isaac test: 103*33e3ecb2SToby Isaac nsize: 4 104*33e3ecb2SToby Isaac suffix: 0 105*33e3ecb2SToby Isaac 106*33e3ecb2SToby Isaac test: 107*33e3ecb2SToby Isaac nsize: 4 108*33e3ecb2SToby Isaac suffix: 0_window 109*33e3ecb2SToby Isaac output_file: output/ex24_0.out 110*33e3ecb2SToby Isaac args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 111*33e3ecb2SToby Isaac requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 112*33e3ecb2SToby Isaac 113*33e3ecb2SToby Isaac TEST*/ 114