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 16d71ae5a4SJacob Faibussowitsch PetscErrorCode PrintVecWithGhosts(DM da, Vec v) 17d71ae5a4SJacob Faibussowitsch { 18c4762a1bSJed Brown PetscScalar **p; 19c4762a1bSJed Brown PetscInt i, j; 20c4762a1bSJed Brown MPI_Comm com; 21c4762a1bSJed Brown PetscMPIInt rank; 22c4762a1bSJed Brown DMDALocalInfo info; 23c4762a1bSJed Brown 24c4762a1bSJed Brown com = PetscObjectComm((PetscObject)da); 259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(com, &rank)); 269566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 2763a3b9bcSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(com, "begin rank %d portion (with ghosts, %" PetscInt_FMT " x %" PetscInt_FMT ")\n", rank, info.gxm, info.gym)); 289566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, v, &p)); 29c4762a1bSJed Brown for (i = info.gxs; i < info.gxs + info.gxm; i++) { 3048a46eb9SPierre Jolivet for (j = info.gys; j < info.gys + info.gym; j++) PetscCall(PetscSynchronizedPrintf(com, "%g, ", (double)PetscRealPart(p[j][i]))); 319566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(com, "\n")); 32c4762a1bSJed Brown } 339566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, v, &p)); 349566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(com, "end rank %d portion\n", rank)); 359566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(com, PETSC_STDOUT)); 363ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 37c4762a1bSJed Brown } 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* Set a Vec v to value, but do not touch ghosts. */ 40d71ae5a4SJacob Faibussowitsch PetscErrorCode VecSetOwned(DM da, Vec v, PetscScalar value) 41d71ae5a4SJacob Faibussowitsch { 42c4762a1bSJed Brown PetscScalar **p; 43c4762a1bSJed Brown PetscInt i, j, xs, xm, ys, ym; 44c4762a1bSJed Brown 459566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, 0, &xm, &ym, 0)); 469566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, v, &p)); 47c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 48ad540459SPierre Jolivet for (j = ys; j < ys + ym; j++) p[j][i] = value; 49c4762a1bSJed Brown } 509566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, v, &p)); 513ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 52c4762a1bSJed Brown } 53c4762a1bSJed Brown 54d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 55d71ae5a4SJacob Faibussowitsch { 56c4762a1bSJed Brown PetscInt M = 4, N = 3; 57c4762a1bSJed Brown DM da; 58c4762a1bSJed Brown Vec local; 59c4762a1bSJed Brown PetscScalar value = 0.0; 60c4762a1bSJed Brown DMBoundaryType bx = DM_BOUNDARY_PERIODIC, by = DM_BOUNDARY_PERIODIC; 61c4762a1bSJed Brown DMDAStencilType stype = DMDA_STENCIL_BOX; 62c4762a1bSJed Brown 63327415f7SBarry Smith PetscFunctionBeginUser; 64*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 659566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD, bx, by, stype, M, N, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da)); 669566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 679566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 689566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(da, &local)); 69c4762a1bSJed Brown 709566063dSJacob Faibussowitsch PetscCall(VecSet(local, value)); 719566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nAfter setting all values to %d:\n", (int)PetscRealPart(value))); 729566063dSJacob Faibussowitsch PetscCall(PrintVecWithGhosts(da, local)); 739566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "done\n")); 74c4762a1bSJed Brown 75c4762a1bSJed Brown value += 1.0; 76c4762a1bSJed Brown /* set values owned by a process, leaving ghosts alone */ 779566063dSJacob Faibussowitsch PetscCall(VecSetOwned(da, local, value)); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* print after re-setting interior values again */ 809566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nAfter setting interior values to %d:\n", (int)PetscRealPart(value))); 819566063dSJacob Faibussowitsch PetscCall(PrintVecWithGhosts(da, local)); 829566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "done\n")); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* communicate ghosts */ 859566063dSJacob Faibussowitsch PetscCall(DMLocalToLocalBegin(da, local, INSERT_VALUES, local)); 869566063dSJacob Faibussowitsch PetscCall(DMLocalToLocalEnd(da, local, INSERT_VALUES, local)); 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* print again */ 899566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nAfter local-to-local communication:\n")); 909566063dSJacob Faibussowitsch PetscCall(PrintVecWithGhosts(da, local)); 919566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "done\n")); 92c4762a1bSJed Brown 93c4762a1bSJed Brown /* Free memory */ 949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&local)); 959566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 969566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 97b122ec5aSJacob Faibussowitsch return 0; 98c4762a1bSJed Brown } 99c4762a1bSJed Brown 100c4762a1bSJed Brown /*TEST 101c4762a1bSJed Brown 102c4762a1bSJed Brown test: 103c4762a1bSJed Brown nsize: 2 104c4762a1bSJed Brown 105c4762a1bSJed Brown TEST*/ 106