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