Actual source code: ex43.c

petsc-3.13.6 2020-09-29
Report Typos and Errors
  1: static char help[] = "Demonstrates the DMLocalToLocal bug in PETSc 3.6.\n\n";

  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

 10:   Contributed by Constantine Khroulev
 11: */

 13:  #include <petscdm.h>
 14:  #include <petscdmda.h>

 16: PetscErrorCode PrintVecWithGhosts(DM da, Vec v)
 17: {
 18:   PetscScalar    **p;
 19:   PetscInt       i, j;
 21:   MPI_Comm       com;
 22:   PetscMPIInt    rank;
 23:   DMDALocalInfo  info;

 25:   com = PetscObjectComm((PetscObject)da);
 26:   MPI_Comm_rank(com, &rank);
 27:   DMDAGetLocalInfo(da, &info);
 28:   PetscSynchronizedPrintf(com, "begin rank %d portion (with ghosts, %D x %D)\n",rank, info.gxm, info.gym);
 29:   DMDAVecGetArray(da, v, &p);
 30:   for (i = info.gxs; i < info.gxs + info.gxm; i++) {
 31:     for (j = info.gys; j < info.gys + info.gym; j++) {
 32:       PetscSynchronizedPrintf(com, "%g, ", (double) PetscRealPart(p[j][i]));
 33:     }
 34:     PetscSynchronizedPrintf(com, "\n");
 35:   }
 36:   DMDAVecRestoreArray(da, v, &p);
 37:   PetscSynchronizedPrintf(com, "end rank %d portion\n", rank);
 38:   PetscSynchronizedFlush(com, PETSC_STDOUT);
 39:   return 0;
 40: }

 42: /* Set a Vec v to value, but do not touch ghosts. */
 43: PetscErrorCode VecSetOwned(DM da, Vec v, PetscScalar value)
 44: {
 45:   PetscScalar    **p;
 46:   PetscInt         i, j, xs, xm, ys, ym;
 47:   PetscErrorCode   ierr;

 49:   DMDAGetCorners(da, &xs, &ys, 0, &xm, &ym, 0);
 50:   DMDAVecGetArray(da, v, &p);
 51:   for (i = xs; i < xs + xm; i++) {
 52:     for (j = ys; j < ys + ym; j++) {
 53:       p[j][i] = value;
 54:     }
 55:   }
 56:   DMDAVecRestoreArray(da, v, &p);
 57:   return 0;
 58: }

 60: int main(int argc, char **argv)
 61: {
 62:   PetscInt         M = 4, N = 3;
 63:   PetscErrorCode   ierr;
 64:   DM               da;
 65:   Vec              local;
 66:   PetscScalar      value = 0.0;
 67:   DMBoundaryType   bx    = DM_BOUNDARY_PERIODIC, by = DM_BOUNDARY_PERIODIC;
 68:   DMDAStencilType  stype = DMDA_STENCIL_BOX;

 70:   PetscInitialize(&argc, &argv, (char*)0, help);if (ierr) return ierr;
 71:   DMDACreate2d(PETSC_COMM_WORLD,bx,by,stype,M,N,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);
 72:   DMSetFromOptions(da);
 73:   DMSetUp(da);
 74:   DMCreateLocalVector(da, &local);

 76:   VecSet(local, value);
 77:   PetscPrintf(PETSC_COMM_WORLD, "\nAfter setting all values to %d:\n", (int)PetscRealPart(value));
 78:   PrintVecWithGhosts(da, local);
 79:   PetscPrintf(PETSC_COMM_WORLD, "done\n");

 81:   value += 1.0;
 82:   /* set values owned by a process, leaving ghosts alone */
 83:   VecSetOwned(da, local, value);

 85:   /* print after re-setting interior values again */
 86:   PetscPrintf(PETSC_COMM_WORLD, "\nAfter setting interior values to %d:\n", (int)PetscRealPart(value));
 87:   PrintVecWithGhosts(da, local);
 88:   PetscPrintf(PETSC_COMM_WORLD, "done\n");

 90:   /* communicate ghosts */
 91:   DMLocalToLocalBegin(da, local, INSERT_VALUES, local);
 92:   DMLocalToLocalEnd(da, local, INSERT_VALUES, local);

 94:   /* print again */
 95:   PetscPrintf(PETSC_COMM_WORLD, "\nAfter local-to-local communication:\n");
 96:   PrintVecWithGhosts(da, local);
 97:   PetscPrintf(PETSC_COMM_WORLD, "done\n");

 99:   /* Free memory */
100:   VecDestroy(&local);
101:   DMDestroy(&da);
102:   PetscFinalize();
103:   return ierr;
104: }


107: /*TEST

109:    test:
110:       nsize: 2

112: TEST*/