Actual source code: ex40.c


  2: static char help[] = "Tests mirror boundary conditions in 2-d.\n\n";

  4: #include <petscdm.h>
  5: #include <petscdmda.h>

  7: int main(int argc, char **argv)
  8: {
  9:   PetscInt       M = 8, N = 8, stencil_width = 1, dof = 1, m, n, xstart, ystart, i, j, c;
 10:   DM             da;
 11:   Vec            global, local;
 12:   PetscScalar ***vglobal;
 13:   PetscViewer    sview;

 15:   PetscFunctionBeginUser;
 16:   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
 17:   PetscCall(PetscOptionsGetInt(NULL, 0, "-stencil_width", &stencil_width, 0));
 18:   PetscCall(PetscOptionsGetInt(NULL, 0, "-dof", &dof, 0));

 20:   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_MIRROR, DM_BOUNDARY_MIRROR, DMDA_STENCIL_STAR, M, N, PETSC_DECIDE, PETSC_DECIDE, dof, stencil_width, NULL, NULL, &da));
 21:   PetscCall(DMSetFromOptions(da));
 22:   PetscCall(DMSetUp(da));
 23:   PetscCall(DMDAGetCorners(da, &xstart, &ystart, 0, &m, &n, 0));

 25:   PetscCall(DMCreateGlobalVector(da, &global));
 26:   PetscCall(DMDAVecGetArrayDOF(da, global, &vglobal));
 27:   for (j = ystart; j < ystart + n; j++) {
 28:     for (i = xstart; i < xstart + m; i++) {
 29:       for (c = 0; c < dof; c++) vglobal[j][i][c] = 100 * j + 10 * (i + 1) + c;
 30:     }
 31:   }
 32:   PetscCall(DMDAVecRestoreArrayDOF(da, global, &vglobal));

 34:   PetscCall(DMCreateLocalVector(da, &local));
 35:   PetscCall(DMGlobalToLocalBegin(da, global, INSERT_VALUES, local));
 36:   PetscCall(DMGlobalToLocalEnd(da, global, INSERT_VALUES, local));

 38:   PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &sview));
 39:   PetscCall(VecView(local, sview));
 40:   PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &sview));
 41:   PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
 42:   PetscCall(VecView(global, PETSC_VIEWER_STDOUT_WORLD));

 44:   PetscCall(DMDestroy(&da));
 45:   PetscCall(VecDestroy(&local));
 46:   PetscCall(VecDestroy(&global));

 48:   PetscCall(PetscFinalize());
 49:   return 0;
 50: }

 52: /*TEST

 54:    test:

 56:    test:
 57:       suffix: 2
 58:       nsize: 4
 59:       filter: grep -v "Vec Object"

 61: TEST*/