Actual source code: ex25.c
2: static char help[] = "Tests DMLocalToGlobal() for dof > 1\n\n";
4: #include <petscdm.h>
5: #include <petscdmda.h>
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;
14: PetscInitialize(&argc,&argv,(char*)0,help);
15: /* Create distributed array and get vectors */
16: 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: DMSetFromOptions(da);
18: DMSetUp(da);
19: DMCreateGlobalVector(da,&global);
20: DMCreateLocalVector(da,&local);
22: DMDAGetCorners(da,&is,&js,&ks,&in,&jen,&kn);
23: 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: DMDAVecRestoreArrayDOF(da,local,&l);
33: DMLocalToGlobalBegin(da,local,ADD_VALUES,global);
34: DMLocalToGlobalEnd(da,local,ADD_VALUES,global);
36: VecView(global,PETSC_VIEWER_STDOUT_WORLD);
38: /* Free memory */
39: VecDestroy(&local);
40: VecDestroy(&global);
41: DMDestroy(&da);
42: PetscFinalize();
43: return 0;
44: }
46: /*TEST
48: test:
49: filter: grep -v -i Process
50: output_file: output/ex25_1.out
52: test:
53: suffix: 2
54: nsize: 2
55: filter: grep -v -i Process
56: output_file: output/ex25_2.out
58: test:
59: suffix: 3
60: nsize: 3
61: filter: grep -v -i Process
62: output_file: output/ex25_2.out
64: test:
65: suffix: 4
66: nsize: 4
67: filter: grep -v -i Process
68: output_file: output/ex25_2.out
70: test:
71: suffix: 5
72: nsize: 5
73: filter: grep -v -i Process
74: output_file: output/ex25_2.out
76: test:
77: suffix: 6
78: nsize: 6
79: filter: grep -v -i Process
80: output_file: output/ex25_2.out
82: TEST*/