Actual source code: ex25.c
petsc-3.7.3 2016-08-01
2: static char help[] = "Tests DMDALocalToGlocal() for dof > 1\n\n";
4: #include <petscdm.h>
5: #include <petscdmda.h>
9: int main(int argc,char **argv)
10: {
11: 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;
13: DM da;
14: Vec local,global;
15: PetscScalar ****l;
17: PetscInitialize(&argc,&argv,(char*)0,help);
19: /* Create distributed array and get vectors */
20: 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);
21: DMCreateGlobalVector(da,&global);
22: DMCreateLocalVector(da,&local);
24: DMDAGetCorners(da,&is,&js,&ks,&in,&jen,&kn);
25: DMDAVecGetArrayDOF(da,local,&l);
26: for (i=is; i<is+in; i++) {
27: for (j=js; j<js+jen; j++) {
28: for (k=ks; k<ks+kn; k++) {
29: l[k][j][i][0] = 2*(i + j*M + k*M*N);
30: l[k][j][i][1] = 2*(i + j*M + k*M*N) + 1;
31: }
32: }
33: }
34: DMDAVecRestoreArrayDOF(da,local,&l);
35: DMLocalToGlobalBegin(da,local,ADD_VALUES,global);
36: DMLocalToGlobalEnd(da,local,ADD_VALUES,global);
38: VecView(global,PETSC_VIEWER_STDOUT_WORLD);
40: /* Free memory */
41: VecDestroy(&local);
42: VecDestroy(&global);
43: DMDestroy(&da);
44: PetscFinalize();
45: return 0;
46: }