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, m = PETSC_DECIDE, n = PETSC_DECIDE, i, j, is, js, in, jen; 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(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, M, N, m, n, 3, 1, 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, 0, &in, &jen, 0)); 24 PetscCall(DMDAVecGetArrayDOF(da, local, &l)); 25 for (i = is; i < is + in; i++) { 26 for (j = js; j < js + jen; j++) { 27 l[j][i][0] = 3 * (i + j * M); 28 l[j][i][1] = 3 * (i + j * M) + 1; 29 l[j][i][2] = 3 * (i + j * M) + 2; 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/ex24_1.out 51 52 test: 53 suffix: 2 54 nsize: 2 55 filter: grep -v -i Process 56 output_file: output/ex24_2.out 57 58 test: 59 suffix: 3 60 nsize: 3 61 filter: grep -v -i Process 62 output_file: output/ex24_2.out 63 64 test: 65 suffix: 4 66 nsize: 4 67 filter: grep -v -i Process 68 output_file: output/ex24_2.out 69 70 test: 71 suffix: 5 72 nsize: 5 73 filter: grep -v -i Process 74 output_file: output/ex24_2.out 75 76 test: 77 suffix: 6 78 nsize: 6 79 filter: grep -v -i Process 80 output_file: output/ex24_2.out 81 82 TEST*/ 83