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++) { 31 PetscCall(PetscSynchronizedPrintf(com, "%g, ", (double) PetscRealPart(p[j][i]))); 32 } 33 PetscCall(PetscSynchronizedPrintf(com, "\n")); 34 } 35 PetscCall(DMDAVecRestoreArray(da, v, &p)); 36 PetscCall(PetscSynchronizedPrintf(com, "end rank %d portion\n", rank)); 37 PetscCall(PetscSynchronizedFlush(com, PETSC_STDOUT)); 38 return 0; 39 } 40 41 /* Set a Vec v to value, but do not touch ghosts. */ 42 PetscErrorCode VecSetOwned(DM da, Vec v, PetscScalar value) 43 { 44 PetscScalar **p; 45 PetscInt i, j, xs, xm, ys, ym; 46 47 PetscCall(DMDAGetCorners(da, &xs, &ys, 0, &xm, &ym, 0)); 48 PetscCall(DMDAVecGetArray(da, v, &p)); 49 for (i = xs; i < xs + xm; i++) { 50 for (j = ys; j < ys + ym; j++) { 51 p[j][i] = value; 52 } 53 } 54 PetscCall(DMDAVecRestoreArray(da, v, &p)); 55 return 0; 56 } 57 58 int main(int argc, char **argv) 59 { 60 PetscInt M = 4, N = 3; 61 DM da; 62 Vec local; 63 PetscScalar value = 0.0; 64 DMBoundaryType bx = DM_BOUNDARY_PERIODIC, by = DM_BOUNDARY_PERIODIC; 65 DMDAStencilType stype = DMDA_STENCIL_BOX; 66 67 PetscCall(PetscInitialize(&argc, &argv, (char*)0, help)); 68 PetscCall(DMDACreate2d(PETSC_COMM_WORLD,bx,by,stype,M,N,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da)); 69 PetscCall(DMSetFromOptions(da)); 70 PetscCall(DMSetUp(da)); 71 PetscCall(DMCreateLocalVector(da, &local)); 72 73 PetscCall(VecSet(local, value)); 74 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nAfter setting all values to %d:\n", (int)PetscRealPart(value))); 75 PetscCall(PrintVecWithGhosts(da, local)); 76 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "done\n")); 77 78 value += 1.0; 79 /* set values owned by a process, leaving ghosts alone */ 80 PetscCall(VecSetOwned(da, local, value)); 81 82 /* print after re-setting interior values again */ 83 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nAfter setting interior values to %d:\n", (int)PetscRealPart(value))); 84 PetscCall(PrintVecWithGhosts(da, local)); 85 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "done\n")); 86 87 /* communicate ghosts */ 88 PetscCall(DMLocalToLocalBegin(da, local, INSERT_VALUES, local)); 89 PetscCall(DMLocalToLocalEnd(da, local, INSERT_VALUES, local)); 90 91 /* print again */ 92 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nAfter local-to-local communication:\n")); 93 PetscCall(PrintVecWithGhosts(da, local)); 94 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "done\n")); 95 96 /* Free memory */ 97 PetscCall(VecDestroy(&local)); 98 PetscCall(DMDestroy(&da)); 99 PetscCall(PetscFinalize()); 100 return 0; 101 } 102 103 /*TEST 104 105 test: 106 nsize: 2 107 108 TEST*/ 109