xref: /petsc/src/dm/tests/ex40.c (revision 66af8762ec03dbef0e079729eb2a1734a35ed7ff)
1 static char help[] = "Tests mirror boundary conditions in 2-d.\n\n";
2 
3 #include <petscdm.h>
4 #include <petscdmda.h>
5 
6 int main(int argc, char **argv)
7 {
8   PetscInt       M = 8, N = 8, stencil_width = 1, dof = 1, m, n, xstart, ystart, i, j, c;
9   DM             da;
10   Vec            global, local;
11   PetscScalar ***vglobal;
12   PetscViewer    sview;
13 
14   PetscFunctionBeginUser;
15   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
16   PetscCall(PetscOptionsGetInt(NULL, 0, "-stencil_width", &stencil_width, 0));
17   PetscCall(PetscOptionsGetInt(NULL, 0, "-dof", &dof, 0));
18 
19   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));
20   PetscCall(DMSetFromOptions(da));
21   PetscCall(DMSetUp(da));
22   PetscCall(DMDAGetCorners(da, &xstart, &ystart, 0, &m, &n, 0));
23 
24   PetscCall(DMCreateGlobalVector(da, &global));
25   PetscCall(DMDAVecGetArrayDOF(da, global, &vglobal));
26   for (j = ystart; j < ystart + n; j++) {
27     for (i = xstart; i < xstart + m; i++) {
28       for (c = 0; c < dof; c++) vglobal[j][i][c] = 100 * j + 10 * (i + 1) + c;
29     }
30   }
31   PetscCall(DMDAVecRestoreArrayDOF(da, global, &vglobal));
32 
33   PetscCall(DMCreateLocalVector(da, &local));
34   PetscCall(DMGlobalToLocalBegin(da, global, INSERT_VALUES, local));
35   PetscCall(DMGlobalToLocalEnd(da, global, INSERT_VALUES, local));
36 
37   PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &sview));
38   PetscCall(VecView(local, sview));
39   PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &sview));
40   PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
41   PetscCall(VecView(global, PETSC_VIEWER_STDOUT_WORLD));
42 
43   PetscCall(DMDestroy(&da));
44   PetscCall(VecDestroy(&local));
45   PetscCall(VecDestroy(&global));
46 
47   PetscCall(PetscFinalize());
48   return 0;
49 }
50 
51 /*TEST
52 
53    test:
54 
55    test:
56       suffix: 2
57       nsize: 4
58       filter: grep -v "Vec Object"
59 
60 TEST*/
61