1 static char help[] = "Demonstrates the DMLocalToLocal bug in PETSc 3.6.\n\n"; 2 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 9 10 Contributed by Constantine Khroulev 11 */ 12 13 #include <petscdm.h> 14 #include <petscdmda.h> 15 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; 23 24 com = PetscObjectComm((PetscObject)da); 25 PetscCallMPI(MPI_Comm_rank(com, &rank)); 26 PetscCall(DMDAGetLocalInfo(da, &info)); 27 PetscCall(PetscSynchronizedPrintf(com, "begin rank %d portion (with ghosts, %" PetscInt_FMT " x %" PetscInt_FMT ")\n", rank, info.gxm, info.gym)); 28 PetscCall(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++) PetscCall(PetscSynchronizedPrintf(com, "%g, ", (double)PetscRealPart(p[j][i]))); 31 PetscCall(PetscSynchronizedPrintf(com, "\n")); 32 } 33 PetscCall(DMDAVecRestoreArray(da, v, &p)); 34 PetscCall(PetscSynchronizedPrintf(com, "end rank %d portion\n", rank)); 35 PetscCall(PetscSynchronizedFlush(com, PETSC_STDOUT)); 36 return PETSC_SUCCESS; 37 } 38 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; 44 45 PetscCall(DMDAGetCorners(da, &xs, &ys, 0, &xm, &ym, 0)); 46 PetscCall(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 PetscCall(DMDAVecRestoreArray(da, v, &p)); 51 return PETSC_SUCCESS; 52 } 53 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; 62 63 PetscFunctionBeginUser; 64 PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 65 PetscCall(DMDACreate2d(PETSC_COMM_WORLD, bx, by, stype, M, N, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da)); 66 PetscCall(DMSetFromOptions(da)); 67 PetscCall(DMSetUp(da)); 68 PetscCall(DMCreateLocalVector(da, &local)); 69 70 PetscCall(VecSet(local, value)); 71 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nAfter setting all values to %d:\n", (int)PetscRealPart(value))); 72 PetscCall(PrintVecWithGhosts(da, local)); 73 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "done\n")); 74 75 value += 1.0; 76 /* set values owned by a process, leaving ghosts alone */ 77 PetscCall(VecSetOwned(da, local, value)); 78 79 /* print after re-setting interior values again */ 80 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nAfter setting interior values to %d:\n", (int)PetscRealPart(value))); 81 PetscCall(PrintVecWithGhosts(da, local)); 82 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "done\n")); 83 84 /* communicate ghosts */ 85 PetscCall(DMLocalToLocalBegin(da, local, INSERT_VALUES, local)); 86 PetscCall(DMLocalToLocalEnd(da, local, INSERT_VALUES, local)); 87 88 /* print again */ 89 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nAfter local-to-local communication:\n")); 90 PetscCall(PrintVecWithGhosts(da, local)); 91 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "done\n")); 92 93 /* Free memory */ 94 PetscCall(VecDestroy(&local)); 95 PetscCall(DMDestroy(&da)); 96 PetscCall(PetscFinalize()); 97 return 0; 98 } 99 100 /*TEST 101 102 test: 103 nsize: 2 104 105 TEST*/ 106