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