1c4762a1bSJed Brown static char help[] = "Demonstrates the DMLocalToLocal bug in PETSc 3.6.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown Use the options 5c4762a1bSJed Brown -da_grid_x <nx> - number of grid points in x direction, if M < 0 6c4762a1bSJed Brown -da_grid_y <ny> - number of grid points in y direction, if N < 0 7c4762a1bSJed Brown -da_processors_x <MX> number of processors in x directio 8c4762a1bSJed Brown -da_processors_y <MY> number of processors in x direction 9c4762a1bSJed Brown 10c4762a1bSJed Brown Contributed by Constantine Khroulev 11c4762a1bSJed Brown */ 12c4762a1bSJed Brown 13c4762a1bSJed Brown #include <petscdm.h> 14c4762a1bSJed Brown #include <petscdmda.h> 15c4762a1bSJed Brown 169371c9d4SSatish Balay PetscErrorCode PrintVecWithGhosts(DM da, Vec v) { 17c4762a1bSJed Brown PetscScalar **p; 18c4762a1bSJed Brown PetscInt i, j; 19c4762a1bSJed Brown MPI_Comm com; 20c4762a1bSJed Brown PetscMPIInt rank; 21c4762a1bSJed Brown DMDALocalInfo info; 22c4762a1bSJed Brown 23c4762a1bSJed Brown com = PetscObjectComm((PetscObject)da); 249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(com, &rank)); 259566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 2663a3b9bcSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(com, "begin rank %d portion (with ghosts, %" PetscInt_FMT " x %" PetscInt_FMT ")\n", rank, info.gxm, info.gym)); 279566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, v, &p)); 28c4762a1bSJed Brown for (i = info.gxs; i < info.gxs + info.gxm; i++) { 29*48a46eb9SPierre Jolivet for (j = info.gys; j < info.gys + info.gym; j++) PetscCall(PetscSynchronizedPrintf(com, "%g, ", (double)PetscRealPart(p[j][i]))); 309566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(com, "\n")); 31c4762a1bSJed Brown } 329566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, v, &p)); 339566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(com, "end rank %d portion\n", rank)); 349566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(com, PETSC_STDOUT)); 35c4762a1bSJed Brown return 0; 36c4762a1bSJed Brown } 37c4762a1bSJed Brown 38c4762a1bSJed Brown /* Set a Vec v to value, but do not touch ghosts. */ 399371c9d4SSatish Balay PetscErrorCode VecSetOwned(DM da, Vec v, PetscScalar value) { 40c4762a1bSJed Brown PetscScalar **p; 41c4762a1bSJed Brown PetscInt i, j, xs, xm, ys, ym; 42c4762a1bSJed Brown 439566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, 0, &xm, &ym, 0)); 449566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, v, &p)); 45c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 469371c9d4SSatish Balay for (j = ys; j < ys + ym; j++) { p[j][i] = value; } 47c4762a1bSJed Brown } 489566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, v, &p)); 49c4762a1bSJed Brown return 0; 50c4762a1bSJed Brown } 51c4762a1bSJed Brown 529371c9d4SSatish Balay int main(int argc, char **argv) { 53c4762a1bSJed Brown PetscInt M = 4, N = 3; 54c4762a1bSJed Brown DM da; 55c4762a1bSJed Brown Vec local; 56c4762a1bSJed Brown PetscScalar value = 0.0; 57c4762a1bSJed Brown DMBoundaryType bx = DM_BOUNDARY_PERIODIC, by = DM_BOUNDARY_PERIODIC; 58c4762a1bSJed Brown DMDAStencilType stype = DMDA_STENCIL_BOX; 59c4762a1bSJed Brown 60327415f7SBarry Smith PetscFunctionBeginUser; 619566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 629566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD, bx, by, stype, M, N, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da)); 639566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 649566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 659566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(da, &local)); 66c4762a1bSJed Brown 679566063dSJacob Faibussowitsch PetscCall(VecSet(local, value)); 689566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nAfter setting all values to %d:\n", (int)PetscRealPart(value))); 699566063dSJacob Faibussowitsch PetscCall(PrintVecWithGhosts(da, local)); 709566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "done\n")); 71c4762a1bSJed Brown 72c4762a1bSJed Brown value += 1.0; 73c4762a1bSJed Brown /* set values owned by a process, leaving ghosts alone */ 749566063dSJacob Faibussowitsch PetscCall(VecSetOwned(da, local, value)); 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* print after re-setting interior values again */ 779566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nAfter setting interior values to %d:\n", (int)PetscRealPart(value))); 789566063dSJacob Faibussowitsch PetscCall(PrintVecWithGhosts(da, local)); 799566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "done\n")); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* communicate ghosts */ 829566063dSJacob Faibussowitsch PetscCall(DMLocalToLocalBegin(da, local, INSERT_VALUES, local)); 839566063dSJacob Faibussowitsch PetscCall(DMLocalToLocalEnd(da, local, INSERT_VALUES, local)); 84c4762a1bSJed Brown 85c4762a1bSJed Brown /* print again */ 869566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nAfter local-to-local communication:\n")); 879566063dSJacob Faibussowitsch PetscCall(PrintVecWithGhosts(da, local)); 889566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "done\n")); 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* Free memory */ 919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&local)); 929566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 939566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 94b122ec5aSJacob Faibussowitsch return 0; 95c4762a1bSJed Brown } 96c4762a1bSJed Brown 97c4762a1bSJed Brown /*TEST 98c4762a1bSJed Brown 99c4762a1bSJed Brown test: 100c4762a1bSJed Brown nsize: 2 101c4762a1bSJed Brown 102c4762a1bSJed Brown TEST*/ 103