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