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