1 static char help[] = "Tests DMLocalToGlobal() for dof > 1\n\n"; 2 3 #include <petscdm.h> 4 #include <petscdmda.h> 5 6 int main(int argc, char **argv) 7 { 8 PetscInt M = 6, N = 5, P = 4, m = PETSC_DECIDE, n = PETSC_DECIDE, p = PETSC_DECIDE, i, j, k, is, js, ks, in, jen, kn; 9 DM da; 10 Vec local, global; 11 PetscScalar ****l; 12 13 PetscFunctionBeginUser; 14 PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 15 /* Create distributed array and get vectors */ 16 PetscCall(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, M, N, P, m, n, p, 2, 1, NULL, NULL, NULL, &da)); 17 PetscCall(DMSetFromOptions(da)); 18 PetscCall(DMSetUp(da)); 19 PetscCall(DMCreateGlobalVector(da, &global)); 20 PetscCall(DMCreateLocalVector(da, &local)); 21 22 PetscCall(DMDAGetCorners(da, &is, &js, &ks, &in, &jen, &kn)); 23 PetscCall(DMDAVecGetArrayDOF(da, local, &l)); 24 for (i = is; i < is + in; i++) { 25 for (j = js; j < js + jen; j++) { 26 for (k = ks; k < ks + kn; k++) { 27 l[k][j][i][0] = 2 * (i + j * M + k * M * N); 28 l[k][j][i][1] = 2 * (i + j * M + k * M * N) + 1; 29 } 30 } 31 } 32 PetscCall(DMDAVecRestoreArrayDOF(da, local, &l)); 33 PetscCall(DMLocalToGlobalBegin(da, local, ADD_VALUES, global)); 34 PetscCall(DMLocalToGlobalEnd(da, local, ADD_VALUES, global)); 35 36 PetscCall(VecView(global, PETSC_VIEWER_STDOUT_WORLD)); 37 38 /* Free memory */ 39 PetscCall(VecDestroy(&local)); 40 PetscCall(VecDestroy(&global)); 41 PetscCall(DMDestroy(&da)); 42 PetscCall(PetscFinalize()); 43 return 0; 44 } 45 46 /*TEST 47 48 test: 49 filter: grep -v -i Process 50 output_file: output/ex25_1.out 51 52 test: 53 suffix: 2 54 nsize: 2 55 filter: grep -v -i Process 56 output_file: output/ex25_2.out 57 58 test: 59 suffix: 3 60 nsize: 3 61 filter: grep -v -i Process 62 output_file: output/ex25_2.out 63 64 test: 65 suffix: 4 66 nsize: 4 67 filter: grep -v -i Process 68 output_file: output/ex25_2.out 69 70 test: 71 suffix: 5 72 nsize: 5 73 filter: grep -v -i Process 74 output_file: output/ex25_2.out 75 76 test: 77 suffix: 6 78 nsize: 6 79 filter: grep -v -i Process 80 output_file: output/ex25_2.out 81 82 TEST*/ 83