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 PetscErrorCode ierr; 11 DM da; 12 Vec local,global; 13 PetscScalar ****l; 14 15 ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 16 /* Create distributed array and get vectors */ 17 ierr = 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);CHKERRQ(ierr); 18 ierr = DMSetFromOptions(da);CHKERRQ(ierr); 19 ierr = DMSetUp(da);CHKERRQ(ierr); 20 ierr = DMCreateGlobalVector(da,&global);CHKERRQ(ierr); 21 ierr = DMCreateLocalVector(da,&local);CHKERRQ(ierr); 22 23 ierr = DMDAGetCorners(da,&is,&js,&ks,&in,&jen,&kn);CHKERRQ(ierr); 24 ierr = DMDAVecGetArrayDOF(da,local,&l);CHKERRQ(ierr); 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 ierr = DMDAVecRestoreArrayDOF(da,local,&l);CHKERRQ(ierr); 34 ierr = DMLocalToGlobalBegin(da,local,ADD_VALUES,global);CHKERRQ(ierr); 35 ierr = DMLocalToGlobalEnd(da,local,ADD_VALUES,global);CHKERRQ(ierr); 36 37 ierr = VecView(global,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 38 39 /* Free memory */ 40 ierr = VecDestroy(&local);CHKERRQ(ierr); 41 ierr = VecDestroy(&global);CHKERRQ(ierr); 42 ierr = DMDestroy(&da);CHKERRQ(ierr); 43 ierr = PetscFinalize(); 44 return ierr; 45 } 46 47 48 49 /*TEST 50 51 test: 52 filter: grep -v -i Process 53 output_file: output/ex25_1.out 54 55 test: 56 suffix: 2 57 nsize: 2 58 filter: grep -v -i Process 59 output_file: output/ex25_2.out 60 61 test: 62 suffix: 3 63 nsize: 3 64 filter: grep -v -i Process 65 output_file: output/ex25_2.out 66 67 test: 68 suffix: 4 69 nsize: 4 70 filter: grep -v -i Process 71 output_file: output/ex25_2.out 72 73 test: 74 suffix: 5 75 nsize: 5 76 filter: grep -v -i Process 77 output_file: output/ex25_2.out 78 79 test: 80 suffix: 6 81 nsize: 6 82 filter: grep -v -i Process 83 output_file: output/ex25_2.out 84 85 TEST*/ 86