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, m = PETSC_DECIDE, n = PETSC_DECIDE, i, j, is, js, in, jen; 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(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, M, N, m, n, 3, 1, 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, 0, &in, &jen, 0)); 23 PetscCall(DMDAVecGetArrayDOF(da, local, &l)); 24 for (i = is; i < is + in; i++) { 25 for (j = js; j < js + jen; j++) { 26 l[j][i][0] = 3 * (i + j * M); 27 l[j][i][1] = 3 * (i + j * M) + 1; 28 l[j][i][2] = 3 * (i + j * M) + 2; 29 } 30 } 31 PetscCall(DMDAVecRestoreArrayDOF(da, local, &l)); 32 PetscCall(DMLocalToGlobalBegin(da, local, ADD_VALUES, global)); 33 PetscCall(DMLocalToGlobalEnd(da, local, ADD_VALUES, global)); 34 35 PetscCall(VecView(global, PETSC_VIEWER_STDOUT_WORLD)); 36 37 /* Free memory */ 38 PetscCall(VecDestroy(&local)); 39 PetscCall(VecDestroy(&global)); 40 PetscCall(DMDestroy(&da)); 41 PetscCall(PetscFinalize()); 42 return 0; 43 } 44 45 /*TEST 46 47 test: 48 filter: grep -v -i Process 49 output_file: output/ex24_1.out 50 51 test: 52 suffix: 2 53 nsize: 2 54 filter: grep -v -i Process 55 output_file: output/ex24_2.out 56 57 test: 58 suffix: 3 59 nsize: 3 60 filter: grep -v -i Process 61 output_file: output/ex24_2.out 62 63 test: 64 suffix: 4 65 nsize: 4 66 filter: grep -v -i Process 67 output_file: output/ex24_2.out 68 69 test: 70 suffix: 5 71 nsize: 5 72 filter: grep -v -i Process 73 output_file: output/ex24_2.out 74 75 test: 76 suffix: 6 77 nsize: 6 78 filter: grep -v -i Process 79 output_file: output/ex24_2.out 80 81 TEST*/ 82