Actual source code: ex4.c

  1: static char help[] = "Test PetscSFFCompose when the ilocal array is not the identity\n\n";

  3: #include <petscsf.h>

  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;

 20:   PetscInitialize(&argc, &argv, NULL, help);
 21:   MPI_Comm_size(PETSC_COMM_WORLD, &size);

 24:   PetscSFCreate(PETSC_COMM_WORLD, &sfA);
 25:   PetscSFCreate(PETSC_COMM_WORLD, &sfB);
 26:   PetscSFSetFromOptions(sfA);
 27:   PetscSFSetFromOptions(sfB);

 29:   PetscOptionsGetBool(NULL, NULL, "-sparse_sfB", &flag, NULL);

 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:   PetscMalloc1(nleavesA, &ilocalA);
 42:   PetscMalloc1(nleavesA, &iremoteA);
 43:   PetscMalloc1(nleavesB, &ilocalB);
 44:   PetscMalloc1(nleavesB, &iremoteB);

 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:   }

 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:   }

 67:   PetscSFSetGraph(sfA, nrootsA, nleavesA, ilocalA, PETSC_OWN_POINTER, iremoteA, PETSC_OWN_POINTER);
 68:   PetscSFSetGraph(sfB, nrootsB, nleavesB, ilocalB, PETSC_OWN_POINTER, iremoteB, PETSC_OWN_POINTER);
 69:   PetscSFSetUp(sfA);
 70:   PetscSFSetUp(sfB);
 71:   PetscObjectSetName((PetscObject)sfA, "sfA");
 72:   PetscObjectSetName((PetscObject)sfB, "sfB");

 74:   VecCreateSeq(PETSC_COMM_WORLD, nrootsA, &a);
 75:   VecCreateSeq(PETSC_COMM_WORLD, nleavesA, &b);
 76:   PetscSFGetLeafRange(sfB, NULL, &maxleafB);
 77:   VecCreateSeq(PETSC_COMM_WORLD, maxleafB + 1, &ba);
 78:   VecGetArray(a, &arrayW);
 79:   for (i = 0; i < nrootsA; i++) arrayW[i] = (PetscScalar)i;
 80:   VecRestoreArray(a, &arrayW);

 82:   PetscPrintf(PETSC_COMM_WORLD, "Initial Vec A\n");
 83:   VecView(a, NULL);
 84:   VecGetArrayRead(a, &arrayR);
 85:   VecGetArray(b, &arrayW);

 87:   PetscSFBcastBegin(sfA, MPIU_SCALAR, arrayR, arrayW, MPI_REPLACE);
 88:   PetscSFBcastEnd(sfA, MPIU_SCALAR, arrayR, arrayW, MPI_REPLACE);
 89:   VecRestoreArray(b, &arrayW);
 90:   VecRestoreArrayRead(a, &arrayR);
 91:   PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast A->B over sfA\n");
 92:   VecView(b, NULL);

 94:   VecGetArrayRead(b, &arrayR);
 95:   VecGetArray(ba, &arrayW);
 96:   arrayW[0] = 10.0; /* Not touched by bcast */
 97:   PetscSFBcastBegin(sfB, MPIU_SCALAR, arrayR, arrayW, MPI_REPLACE);
 98:   PetscSFBcastEnd(sfB, MPIU_SCALAR, arrayR, arrayW, MPI_REPLACE);
 99:   VecRestoreArrayRead(b, &arrayR);
100:   VecRestoreArray(ba, &arrayW);

102:   PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast B->BA over sfB\n");
103:   VecView(ba, NULL);

105:   PetscSFCompose(sfA, sfB, &sfBA);
106:   PetscSFSetFromOptions(sfBA);
107:   PetscObjectSetName((PetscObject)sfBA, "(sfB o sfA)");
108:   VecGetArrayRead(a, &arrayR);
109:   VecGetArray(ba, &arrayW);
110:   arrayW[0] = 11.0; /* Not touched by bcast */
111:   PetscSFBcastBegin(sfBA, MPIU_SCALAR, arrayR, arrayW, MPI_REPLACE);
112:   PetscSFBcastEnd(sfBA, MPIU_SCALAR, arrayR, arrayW, MPI_REPLACE);
113:   VecRestoreArray(ba, &arrayW);
114:   VecRestoreArrayRead(a, &arrayR);
115:   PetscPrintf(PETSC_COMM_WORLD, "\nBroadcast A->BA over sfBA (sfB o sfA)\n");
116:   VecView(ba, NULL);

118:   VecDestroy(&ba);
119:   VecDestroy(&b);
120:   VecDestroy(&a);

122:   PetscSFView(sfA, NULL);
123:   PetscSFView(sfB, NULL);
124:   PetscSFView(sfBA, NULL);
125:   PetscSFDestroy(&sfA);
126:   PetscSFDestroy(&sfB);
127:   PetscSFDestroy(&sfBA);

129:   PetscFinalize();
130:   return 0;
131: }

133: /*TEST

135:    test:
136:      suffix: 1

138:    test:
139:      suffix: 2
140:      filter: grep -v "type" | grep -v "sort"
141:      args: -sparse_sfB

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)

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)

158: TEST*/