Actual source code: ex43.c

  1: static char help[] = "Demonstrates the DMLocalToLocal bug in PETSc 3.6.\n\n";

  3: /*
  4: Use the options
  5:      -da_grid_x <nx> - number of grid points in x direction, if M < 0
  6:      -da_grid_y <ny> - number of grid points in y direction, if N < 0
  7:      -da_processors_x <MX> number of processors in x directio
  8:      -da_processors_y <MY> number of processors in x direction

 10:   Contributed by Constantine Khroulev
 11: */

 13: #include <petscdm.h>
 14: #include <petscdmda.h>

 16: PetscErrorCode PrintVecWithGhosts(DM da, Vec v)
 17: {
 18:   PetscScalar **p;
 19:   PetscInt      i, j;
 20:   MPI_Comm      com;
 21:   PetscMPIInt   rank;
 22:   DMDALocalInfo info;

 24:   com = PetscObjectComm((PetscObject)da);
 25:   MPI_Comm_rank(com, &rank);
 26:   DMDAGetLocalInfo(da, &info);
 27:   PetscSynchronizedPrintf(com, "begin rank %d portion (with ghosts, %" PetscInt_FMT " x %" PetscInt_FMT ")\n", rank, info.gxm, info.gym);
 28:   DMDAVecGetArray(da, v, &p);
 29:   for (i = info.gxs; i < info.gxs + info.gxm; i++) {
 30:     for (j = info.gys; j < info.gys + info.gym; j++) PetscSynchronizedPrintf(com, "%g, ", (double)PetscRealPart(p[j][i]));
 31:     PetscSynchronizedPrintf(com, "\n");
 32:   }
 33:   DMDAVecRestoreArray(da, v, &p);
 34:   PetscSynchronizedPrintf(com, "end rank %d portion\n", rank);
 35:   PetscSynchronizedFlush(com, PETSC_STDOUT);
 36:   return 0;
 37: }

 39: /* Set a Vec v to value, but do not touch ghosts. */
 40: PetscErrorCode VecSetOwned(DM da, Vec v, PetscScalar value)
 41: {
 42:   PetscScalar **p;
 43:   PetscInt      i, j, xs, xm, ys, ym;

 45:   DMDAGetCorners(da, &xs, &ys, 0, &xm, &ym, 0);
 46:   DMDAVecGetArray(da, v, &p);
 47:   for (i = xs; i < xs + xm; i++) {
 48:     for (j = ys; j < ys + ym; j++) p[j][i] = value;
 49:   }
 50:   DMDAVecRestoreArray(da, v, &p);
 51:   return 0;
 52: }

 54: int main(int argc, char **argv)
 55: {
 56:   PetscInt        M = 4, N = 3;
 57:   DM              da;
 58:   Vec             local;
 59:   PetscScalar     value = 0.0;
 60:   DMBoundaryType  bx = DM_BOUNDARY_PERIODIC, by = DM_BOUNDARY_PERIODIC;
 61:   DMDAStencilType stype = DMDA_STENCIL_BOX;

 64:   PetscInitialize(&argc, &argv, (char *)0, help);
 65:   DMDACreate2d(PETSC_COMM_WORLD, bx, by, stype, M, N, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da);
 66:   DMSetFromOptions(da);
 67:   DMSetUp(da);
 68:   DMCreateLocalVector(da, &local);

 70:   VecSet(local, value);
 71:   PetscPrintf(PETSC_COMM_WORLD, "\nAfter setting all values to %d:\n", (int)PetscRealPart(value));
 72:   PrintVecWithGhosts(da, local);
 73:   PetscPrintf(PETSC_COMM_WORLD, "done\n");

 75:   value += 1.0;
 76:   /* set values owned by a process, leaving ghosts alone */
 77:   VecSetOwned(da, local, value);

 79:   /* print after re-setting interior values again */
 80:   PetscPrintf(PETSC_COMM_WORLD, "\nAfter setting interior values to %d:\n", (int)PetscRealPart(value));
 81:   PrintVecWithGhosts(da, local);
 82:   PetscPrintf(PETSC_COMM_WORLD, "done\n");

 84:   /* communicate ghosts */
 85:   DMLocalToLocalBegin(da, local, INSERT_VALUES, local);
 86:   DMLocalToLocalEnd(da, local, INSERT_VALUES, local);

 88:   /* print again */
 89:   PetscPrintf(PETSC_COMM_WORLD, "\nAfter local-to-local communication:\n");
 90:   PrintVecWithGhosts(da, local);
 91:   PetscPrintf(PETSC_COMM_WORLD, "done\n");

 93:   /* Free memory */
 94:   VecDestroy(&local);
 95:   DMDestroy(&da);
 96:   PetscFinalize();
 97:   return 0;
 98: }

100: /*TEST

102:    test:
103:       nsize: 2

105: TEST*/