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