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