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 16c4762a1bSJed Brown PetscErrorCode PrintVecWithGhosts(DM da, Vec v) 17c4762a1bSJed Brown { 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); 255f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(com, &rank)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(da, &info)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(com, "begin rank %d portion (with ghosts, %D x %D)\n",rank, info.gxm, info.gym)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da, v, &p)); 29c4762a1bSJed Brown for (i = info.gxs; i < info.gxs + info.gxm; i++) { 30c4762a1bSJed Brown for (j = info.gys; j < info.gys + info.gym; j++) { 315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(com, "%g, ", (double) PetscRealPart(p[j][i]))); 32c4762a1bSJed Brown } 335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(com, "\n")); 34c4762a1bSJed Brown } 355f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da, v, &p)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(com, "end rank %d portion\n", rank)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(com, PETSC_STDOUT)); 38c4762a1bSJed Brown return 0; 39c4762a1bSJed Brown } 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* Set a Vec v to value, but do not touch ghosts. */ 42c4762a1bSJed Brown PetscErrorCode VecSetOwned(DM da, Vec v, PetscScalar value) 43c4762a1bSJed Brown { 44c4762a1bSJed Brown PetscScalar **p; 45c4762a1bSJed Brown PetscInt i, j, xs, xm, ys, ym; 46c4762a1bSJed Brown 475f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da, &xs, &ys, 0, &xm, &ym, 0)); 485f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da, v, &p)); 49c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 50c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 51c4762a1bSJed Brown p[j][i] = value; 52c4762a1bSJed Brown } 53c4762a1bSJed Brown } 545f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da, v, &p)); 55c4762a1bSJed Brown return 0; 56c4762a1bSJed Brown } 57c4762a1bSJed Brown 58c4762a1bSJed Brown int main(int argc, char **argv) 59c4762a1bSJed Brown { 60c4762a1bSJed Brown PetscInt M = 4, N = 3; 61c4762a1bSJed Brown DM da; 62c4762a1bSJed Brown Vec local; 63c4762a1bSJed Brown PetscScalar value = 0.0; 64c4762a1bSJed Brown DMBoundaryType bx = DM_BOUNDARY_PERIODIC, by = DM_BOUNDARY_PERIODIC; 65c4762a1bSJed Brown DMDAStencilType stype = DMDA_STENCIL_BOX; 66c4762a1bSJed Brown 67*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, (char*)0, help)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,bx,by,stype,M,N,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLocalVector(da, &local)); 72c4762a1bSJed Brown 735f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(local, value)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\nAfter setting all values to %d:\n", (int)PetscRealPart(value))); 755f80ce2aSJacob Faibussowitsch CHKERRQ(PrintVecWithGhosts(da, local)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "done\n")); 77c4762a1bSJed Brown 78c4762a1bSJed Brown value += 1.0; 79c4762a1bSJed Brown /* set values owned by a process, leaving ghosts alone */ 805f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetOwned(da, local, value)); 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* print after re-setting interior values again */ 835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\nAfter setting interior values to %d:\n", (int)PetscRealPart(value))); 845f80ce2aSJacob Faibussowitsch CHKERRQ(PrintVecWithGhosts(da, local)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "done\n")); 86c4762a1bSJed Brown 87c4762a1bSJed Brown /* communicate ghosts */ 885f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToLocalBegin(da, local, INSERT_VALUES, local)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToLocalEnd(da, local, INSERT_VALUES, local)); 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* print again */ 925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\nAfter local-to-local communication:\n")); 935f80ce2aSJacob Faibussowitsch CHKERRQ(PrintVecWithGhosts(da, local)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "done\n")); 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* Free memory */ 975f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&local)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 99*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 100*b122ec5aSJacob Faibussowitsch return 0; 101c4762a1bSJed Brown } 102c4762a1bSJed Brown 103c4762a1bSJed Brown /*TEST 104c4762a1bSJed Brown 105c4762a1bSJed Brown test: 106c4762a1bSJed Brown nsize: 2 107c4762a1bSJed Brown 108c4762a1bSJed Brown TEST*/ 109