Actual source code: ex41.c
2: static char help[] = "Tests mirror boundary conditions in 3-d.\n\n";
4: #include <petscdm.h>
5: #include <petscdmda.h>
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;
17: PetscInitialize(&argc, &argv, (char *)0, help);
18: PetscOptionsGetInt(NULL, 0, "-stencil_width", &stencil_width, 0);
19: PetscOptionsGetInt(NULL, 0, "-dof", &dof, 0);
21: 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: DMSetFromOptions(da);
23: DMSetUp(da);
24: DMDAGetCorners(da, &xstart, &ystart, &zstart, &m, &n, &p);
26: DMCreateGlobalVector(da, &global);
27: 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: DMDAVecRestoreArrayDOF(da, global, &vglobal);
37: DMCreateLocalVector(da, &local);
38: DMGlobalToLocalBegin(da, global, ADD_VALUES, local);
39: DMGlobalToLocalEnd(da, global, ADD_VALUES, local);
41: VecSum(local, &sum);
42: PetscSynchronizedPrintf(PETSC_COMM_WORLD, "sum %g\n", (double)PetscRealPart(sum));
43: PetscSynchronizedFlush(PETSC_COMM_WORLD, stdout);
44: PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &sview);
45: VecView(local, sview);
46: PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &sview);
47: PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);
49: DMDestroy(&da);
50: VecDestroy(&local);
51: VecDestroy(&global);
53: PetscFinalize();
54: return 0;
55: }
57: /*TEST
59: test:
61: test:
62: suffix: 2
63: nsize: 3
65: TEST*/