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