1 static char help[] = "Test PetscSFFCompose when the ilocal array is not the identity\n\n"; 2 3 #include <petscsf.h> 4 5 int main(int argc, char **argv) 6 { 7 PetscSF sfA, sfB, sfBA; 8 PetscInt nrootsA, nleavesA, nrootsB, nleavesB; 9 PetscInt *ilocalA, *ilocalB; 10 PetscSFNode *iremoteA, *iremoteB; 11 Vec a, b, ba; 12 const PetscScalar *arrayR; 13 PetscScalar *arrayW; 14 PetscMPIInt size; 15 PetscInt i; 16 PetscInt maxleafB; 17 PetscBool flag = PETSC_FALSE; 18 19 PetscFunctionBeginUser; 20 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 21 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 22 PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_USER, "Only coded for one MPI process"); 23 24 PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sfA)); 25 PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sfB)); 26 PetscCall(PetscSFSetFromOptions(sfA)); 27 PetscCall(PetscSFSetFromOptions(sfB)); 28 29 PetscCall(PetscOptionsGetBool(NULL, NULL, "-sparse_sfB", &flag, NULL)); 30 31 if (flag) { 32 /* sfA permutes indices, sfB has sparse leaf space. */ 33 nrootsA = 3; 34 nleavesA = 3; 35 nrootsB = 3; 36 nleavesB = 2; 37 } else { 38 /* sfA reverses indices, sfB is identity */ 39 nrootsA = nrootsB = nleavesA = nleavesB = 4; 40 } 41 PetscCall(PetscMalloc1(nleavesA, &ilocalA)); 42 PetscCall(PetscMalloc1(nleavesA, &iremoteA)); 43 PetscCall(PetscMalloc1(nleavesB, &ilocalB)); 44 PetscCall(PetscMalloc1(nleavesB, &iremoteB)); 45 46 for (i = 0; i < nleavesA; i++) { 47 iremoteA[i].rank = 0; 48 iremoteA[i].index = i; 49 if (flag) { 50 ilocalA[i] = (i + 1) % nleavesA; 51 } else { 52 ilocalA[i] = nleavesA - i - 1; 53 } 54 } 55 56 for (i = 0; i < nleavesB; i++) { 57 iremoteB[i].rank = 0; 58 if (flag) { 59 ilocalB[i] = nleavesB - i; 60 iremoteB[i].index = nleavesB - i - 1; 61 } else { 62 ilocalB[i] = i; 63 iremoteB[i].index = i; 64 } 65 } 66 67 PetscCall(PetscSFSetGraph(sfA, nrootsA, nleavesA, ilocalA, PETSC_OWN_POINTER, iremoteA, PETSC_OWN_POINTER)); 68 PetscCall(PetscSFSetGraph(sfB, nrootsB, nleavesB, ilocalB, PETSC_OWN_POINTER, iremoteB, PETSC_OWN_POINTER)); 69 PetscCall(PetscSFSetUp(sfA)); 70 PetscCall(PetscSFSetUp(sfB)); 71 PetscCall(PetscObjectSetName((PetscObject)sfA, "sfA")); 72 PetscCall(PetscObjectSetName((PetscObject)sfB, "sfB")); 73 74 PetscCall(VecCreateSeq(PETSC_COMM_WORLD, nrootsA, &a)); 75 PetscCall(VecCreateSeq(PETSC_COMM_WORLD, nleavesA, &b)); 76 PetscCall(PetscSFGetLeafRange(sfB, NULL, &maxleafB)); 77 PetscCall(VecCreateSeq(PETSC_COMM_WORLD, maxleafB + 1, &ba)); 78 PetscCall(VecGetArray(a, &arrayW)); 79 for (i = 0; i < nrootsA; i++) arrayW[i] = (PetscScalar)i; 80 PetscCall(VecRestoreArray(a, &arrayW)); 81 82 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial Vec A\n")); 83 PetscCall(VecView(a, NULL)); 84 PetscCall(VecGetArrayRead(a, &arrayR)); 85 PetscCall(VecGetArray(b, &arrayW)); 86 87 PetscCall(PetscSFBcastBegin(sfA, MPIU_SCALAR, arrayR, arrayW, MPI_REPLACE)); 88 PetscCall(PetscSFBcastEnd(sfA, MPIU_SCALAR, arrayR, arrayW, MPI_REPLACE)); 89 PetscCall(VecRestoreArray(b, &arrayW)); 90 PetscCall(VecRestoreArrayRead(a, &arrayR)); 91 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast A->B over sfA\n")); 92 PetscCall(VecView(b, NULL)); 93 94 PetscCall(VecGetArrayRead(b, &arrayR)); 95 PetscCall(VecGetArray(ba, &arrayW)); 96 arrayW[0] = 10.0; /* Not touched by bcast */ 97 PetscCall(PetscSFBcastBegin(sfB, MPIU_SCALAR, arrayR, arrayW, MPI_REPLACE)); 98 PetscCall(PetscSFBcastEnd(sfB, MPIU_SCALAR, arrayR, arrayW, MPI_REPLACE)); 99 PetscCall(VecRestoreArrayRead(b, &arrayR)); 100 PetscCall(VecRestoreArray(ba, &arrayW)); 101 102 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast B->BA over sfB\n")); 103 PetscCall(VecView(ba, NULL)); 104 105 PetscCall(PetscSFCompose(sfA, sfB, &sfBA)); 106 PetscCall(PetscSFSetFromOptions(sfBA)); 107 PetscCall(PetscObjectSetName((PetscObject)sfBA, "(sfB o sfA)")); 108 PetscCall(VecGetArrayRead(a, &arrayR)); 109 PetscCall(VecGetArray(ba, &arrayW)); 110 arrayW[0] = 11.0; /* Not touched by bcast */ 111 PetscCall(PetscSFBcastBegin(sfBA, MPIU_SCALAR, arrayR, arrayW, MPI_REPLACE)); 112 PetscCall(PetscSFBcastEnd(sfBA, MPIU_SCALAR, arrayR, arrayW, MPI_REPLACE)); 113 PetscCall(VecRestoreArray(ba, &arrayW)); 114 PetscCall(VecRestoreArrayRead(a, &arrayR)); 115 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast A->BA over sfBA (sfB o sfA)\n")); 116 PetscCall(VecView(ba, NULL)); 117 118 PetscCall(VecDestroy(&ba)); 119 PetscCall(VecDestroy(&b)); 120 PetscCall(VecDestroy(&a)); 121 122 PetscCall(PetscSFView(sfA, NULL)); 123 PetscCall(PetscSFView(sfB, NULL)); 124 PetscCall(PetscSFView(sfBA, NULL)); 125 PetscCall(PetscSFDestroy(&sfA)); 126 PetscCall(PetscSFDestroy(&sfB)); 127 PetscCall(PetscSFDestroy(&sfBA)); 128 129 PetscCall(PetscFinalize()); 130 return 0; 131 } 132 133 /*TEST 134 135 test: 136 suffix: 1 137 138 test: 139 suffix: 2 140 filter: grep -v "type" | grep -v "sort" 141 args: -sparse_sfB 142 143 test: 144 suffix: 2_window 145 filter: grep -v "type" | grep -v "sort" 146 output_file: output/ex4_2.out 147 args: -sparse_sfB -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 148 requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 149 150 # The nightly test suite with MPICH uses ch3:sock, which is broken when winsize == 0 in some of the processes 151 test: 152 suffix: 2_window_shared 153 filter: grep -v "type" | grep -v "sort" 154 output_file: output/ex4_2.out 155 args: -sparse_sfB -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor shared 156 requires: defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) !defined(PETSC_HAVE_MPICH_NUMVERSION) defined(PETSC_HAVE_MPI_ONE_SIDED) 157 158 TEST*/ 159