xref: /petsc/src/dm/tests/ex43.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
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   CHKERRMPI(MPI_Comm_rank(com, &rank));
26   CHKERRQ(DMDAGetLocalInfo(da, &info));
27   CHKERRQ(PetscSynchronizedPrintf(com, "begin rank %d portion (with ghosts, %D x %D)\n",rank, info.gxm, info.gym));
28   CHKERRQ(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       CHKERRQ(PetscSynchronizedPrintf(com, "%g, ", (double) PetscRealPart(p[j][i])));
32     }
33     CHKERRQ(PetscSynchronizedPrintf(com, "\n"));
34   }
35   CHKERRQ(DMDAVecRestoreArray(da, v, &p));
36   CHKERRQ(PetscSynchronizedPrintf(com, "end rank %d portion\n", rank));
37   CHKERRQ(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   CHKERRQ(DMDAGetCorners(da, &xs, &ys, 0, &xm, &ym, 0));
48   CHKERRQ(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   CHKERRQ(DMDAVecRestoreArray(da, v, &p));
55   return 0;
56 }
57 
58 int main(int argc, char **argv)
59 {
60   PetscInt         M = 4, N = 3;
61   PetscErrorCode   ierr;
62   DM               da;
63   Vec              local;
64   PetscScalar      value = 0.0;
65   DMBoundaryType   bx    = DM_BOUNDARY_PERIODIC, by = DM_BOUNDARY_PERIODIC;
66   DMDAStencilType  stype = DMDA_STENCIL_BOX;
67 
68   ierr = PetscInitialize(&argc, &argv, (char*)0, help);if (ierr) return ierr;
69   CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,bx,by,stype,M,N,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da));
70   CHKERRQ(DMSetFromOptions(da));
71   CHKERRQ(DMSetUp(da));
72   CHKERRQ(DMCreateLocalVector(da, &local));
73 
74   CHKERRQ(VecSet(local, value));
75   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\nAfter setting all values to %d:\n", (int)PetscRealPart(value)));
76   CHKERRQ(PrintVecWithGhosts(da, local));
77   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "done\n"));
78 
79   value += 1.0;
80   /* set values owned by a process, leaving ghosts alone */
81   CHKERRQ(VecSetOwned(da, local, value));
82 
83   /* print after re-setting interior values again */
84   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\nAfter setting interior values to %d:\n", (int)PetscRealPart(value)));
85   CHKERRQ(PrintVecWithGhosts(da, local));
86   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "done\n"));
87 
88   /* communicate ghosts */
89   CHKERRQ(DMLocalToLocalBegin(da, local, INSERT_VALUES, local));
90   CHKERRQ(DMLocalToLocalEnd(da, local, INSERT_VALUES, local));
91 
92   /* print again */
93   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\nAfter local-to-local communication:\n"));
94   CHKERRQ(PrintVecWithGhosts(da, local));
95   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "done\n"));
96 
97   /* Free memory */
98   CHKERRQ(VecDestroy(&local));
99   CHKERRQ(DMDestroy(&da));
100   ierr = PetscFinalize();
101   return ierr;
102 }
103 
104 /*TEST
105 
106    test:
107       nsize: 2
108 
109 TEST*/
110