xref: /petsc/src/dm/tests/ex43.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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);
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++) {
30c4762a1bSJed Brown     for (j = info.gys; j < info.gys + info.gym; j++) {
319566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedPrintf(com, "%g, ", (double) PetscRealPart(p[j][i])));
32c4762a1bSJed Brown     }
339566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(com, "\n"));
34c4762a1bSJed Brown   }
359566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, v, &p));
369566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedPrintf(com, "end rank %d portion\n", rank));
379566063dSJacob Faibussowitsch   PetscCall(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 
479566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, 0, &xm, &ym, 0));
489566063dSJacob Faibussowitsch   PetscCall(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   }
549566063dSJacob Faibussowitsch   PetscCall(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*327415f7SBarry Smith   PetscFunctionBeginUser;
689566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char*)0, help));
699566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(PETSC_COMM_WORLD,bx,by,stype,M,N,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da));
709566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
719566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
729566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(da, &local));
73c4762a1bSJed Brown 
749566063dSJacob Faibussowitsch   PetscCall(VecSet(local, value));
759566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nAfter setting all values to %d:\n", (int)PetscRealPart(value)));
769566063dSJacob Faibussowitsch   PetscCall(PrintVecWithGhosts(da, local));
779566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "done\n"));
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   value += 1.0;
80c4762a1bSJed Brown   /* set values owned by a process, leaving ghosts alone */
819566063dSJacob Faibussowitsch   PetscCall(VecSetOwned(da, local, value));
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   /* print after re-setting interior values again */
849566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nAfter setting interior values to %d:\n", (int)PetscRealPart(value)));
859566063dSJacob Faibussowitsch   PetscCall(PrintVecWithGhosts(da, local));
869566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "done\n"));
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   /* communicate ghosts */
899566063dSJacob Faibussowitsch   PetscCall(DMLocalToLocalBegin(da, local, INSERT_VALUES, local));
909566063dSJacob Faibussowitsch   PetscCall(DMLocalToLocalEnd(da, local, INSERT_VALUES, local));
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   /* print again */
939566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nAfter local-to-local communication:\n"));
949566063dSJacob Faibussowitsch   PetscCall(PrintVecWithGhosts(da, local));
959566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "done\n"));
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   /* Free memory */
989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&local));
999566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1009566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
101b122ec5aSJacob Faibussowitsch   return 0;
102c4762a1bSJed Brown }
103c4762a1bSJed Brown 
104c4762a1bSJed Brown /*TEST
105c4762a1bSJed Brown 
106c4762a1bSJed Brown    test:
107c4762a1bSJed Brown       nsize: 2
108c4762a1bSJed Brown 
109c4762a1bSJed Brown TEST*/
110