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 PetscErrorCode ierr; 8 PetscSF sfA, sfB, sfBA; 9 PetscInt nrootsA, nleavesA, nrootsB, nleavesB; 10 PetscInt *ilocalA, *ilocalB; 11 PetscSFNode *iremoteA, *iremoteB; 12 Vec a, b, ba; 13 const PetscScalar *arrayR; 14 PetscScalar *arrayW; 15 PetscMPIInt size; 16 PetscInt i; 17 PetscInt maxleafB; 18 PetscBool flag = PETSC_FALSE; 19 20 ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 21 CHKERRMPI(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 CHKERRQ(PetscSFCreate(PETSC_COMM_WORLD, &sfA)); 25 CHKERRQ(PetscSFCreate(PETSC_COMM_WORLD, &sfB)); 26 CHKERRQ(PetscSFSetFromOptions(sfA)); 27 CHKERRQ(PetscSFSetFromOptions(sfB)); 28 29 CHKERRQ(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 CHKERRQ(PetscMalloc1(nleavesA, &ilocalA)); 42 CHKERRQ(PetscMalloc1(nleavesA, &iremoteA)); 43 CHKERRQ(PetscMalloc1(nleavesB, &ilocalB)); 44 CHKERRQ(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 CHKERRQ(PetscSFSetGraph(sfA, nrootsA, nleavesA, ilocalA, PETSC_OWN_POINTER, iremoteA, PETSC_OWN_POINTER)); 68 CHKERRQ(PetscSFSetGraph(sfB, nrootsB, nleavesB, ilocalB, PETSC_OWN_POINTER, iremoteB, PETSC_OWN_POINTER)); 69 CHKERRQ(PetscSFSetUp(sfA)); 70 CHKERRQ(PetscSFSetUp(sfB)); 71 CHKERRQ(PetscObjectSetName((PetscObject)sfA, "sfA")); 72 CHKERRQ(PetscObjectSetName((PetscObject)sfB, "sfB")); 73 74 CHKERRQ(VecCreateSeq(PETSC_COMM_WORLD, nrootsA, &a)); 75 CHKERRQ(VecCreateSeq(PETSC_COMM_WORLD, nleavesA, &b)); 76 CHKERRQ(PetscSFGetLeafRange(sfB, NULL, &maxleafB)); 77 CHKERRQ(VecCreateSeq(PETSC_COMM_WORLD, maxleafB+1, &ba)); 78 CHKERRQ(VecGetArray(a, &arrayW)); 79 for (i = 0; i < nrootsA; i++) { 80 arrayW[i] = (PetscScalar)i; 81 } 82 CHKERRQ(VecRestoreArray(a, &arrayW)); 83 84 CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Initial Vec A\n")); 85 CHKERRQ(VecView(a, NULL)); 86 CHKERRQ(VecGetArrayRead(a, &arrayR)); 87 CHKERRQ(VecGetArray(b, &arrayW)); 88 89 CHKERRQ(PetscSFBcastBegin(sfA, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE)); 90 CHKERRQ(PetscSFBcastEnd(sfA, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE)); 91 CHKERRQ(VecRestoreArray(b, &arrayW)); 92 CHKERRQ(VecRestoreArrayRead(a, &arrayR)); 93 CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast A->B over sfA\n")); 94 CHKERRQ(VecView(b, NULL)); 95 96 CHKERRQ(VecGetArrayRead(b, &arrayR)); 97 CHKERRQ(VecGetArray(ba, &arrayW)); 98 arrayW[0] = 10.0; /* Not touched by bcast */ 99 CHKERRQ(PetscSFBcastBegin(sfB, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE)); 100 CHKERRQ(PetscSFBcastEnd(sfB, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE)); 101 CHKERRQ(VecRestoreArrayRead(b, &arrayR)); 102 CHKERRQ(VecRestoreArray(ba, &arrayW)); 103 104 CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast B->BA over sfB\n")); 105 CHKERRQ(VecView(ba, NULL)); 106 107 CHKERRQ(PetscSFCompose(sfA, sfB, &sfBA)); 108 CHKERRQ(PetscSFSetFromOptions(sfBA)); 109 CHKERRQ(PetscObjectSetName((PetscObject)sfBA, "(sfB o sfA)")); 110 CHKERRQ(VecGetArrayRead(a, &arrayR)); 111 CHKERRQ(VecGetArray(ba, &arrayW)); 112 arrayW[0] = 11.0; /* Not touched by bcast */ 113 CHKERRQ(PetscSFBcastBegin(sfBA, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE)); 114 CHKERRQ(PetscSFBcastEnd(sfBA, MPIU_SCALAR, arrayR, arrayW,MPI_REPLACE)); 115 CHKERRQ(VecRestoreArray(ba, &arrayW)); 116 CHKERRQ(VecRestoreArrayRead(a, &arrayR)); 117 CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast A->BA over sfBA (sfB o sfA)\n")); 118 CHKERRQ(VecView(ba, NULL)); 119 120 CHKERRQ(VecDestroy(&ba)); 121 CHKERRQ(VecDestroy(&b)); 122 CHKERRQ(VecDestroy(&a)); 123 124 CHKERRQ(PetscSFView(sfA, NULL)); 125 CHKERRQ(PetscSFView(sfB, NULL)); 126 CHKERRQ(PetscSFView(sfBA, NULL)); 127 CHKERRQ(PetscSFDestroy(&sfA)); 128 CHKERRQ(PetscSFDestroy(&sfB)); 129 CHKERRQ(PetscSFDestroy(&sfBA)); 130 131 ierr = PetscFinalize(); 132 return ierr; 133 } 134 135 /*TEST 136 137 test: 138 suffix: 1 139 140 test: 141 suffix: 2 142 filter: grep -v "type" | grep -v "sort" 143 args: -sparse_sfB 144 145 test: 146 suffix: 2_window 147 filter: grep -v "type" | grep -v "sort" 148 output_file: output/ex4_2.out 149 args: -sparse_sfB -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}} 150 requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 151 152 # The nightly test suite with MPICH uses ch3:sock, which is broken when winsize == 0 in some of the processes 153 test: 154 suffix: 2_window_shared 155 filter: grep -v "type" | grep -v "sort" 156 output_file: output/ex4_2.out 157 args: -sparse_sfB -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor shared 158 requires: defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) !defined(PETSC_HAVE_MPICH_NUMVERSION) defined(PETSC_HAVE_MPI_ONE_SIDED) 159 160 TEST*/ 161