1 
2 static char help[] = "Tests mirror boundary conditions in 2-d.\n\n";
3 
4 #include <petscdm.h>
5 #include <petscdmda.h>
6 
main(int argc,char ** argv)7 int main(int argc,char **argv)
8 {
9   PetscErrorCode ierr;
10   PetscInt       M = 8, N = 8,stencil_width = 1, dof = 1,m,n,xstart,ystart,i,j,c;
11   DM             da;
12   Vec            global,local;
13   PetscScalar    ***vglobal;
14   PetscViewer    sview;
15 
16   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
17   ierr = PetscOptionsGetInt(NULL,0,"-stencil_width",&stencil_width,0);CHKERRQ(ierr);
18   ierr = PetscOptionsGetInt(NULL,0,"-dof",&dof,0);CHKERRQ(ierr);
19 
20   ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_MIRROR,DM_BOUNDARY_MIRROR,DMDA_STENCIL_STAR,M,N,PETSC_DECIDE,PETSC_DECIDE,dof,stencil_width,NULL,NULL,&da);CHKERRQ(ierr);
21   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
22   ierr = DMSetUp(da);CHKERRQ(ierr);
23   ierr = DMDAGetCorners(da,&xstart,&ystart,0,&m,&n,0);CHKERRQ(ierr);
24 
25   ierr = DMCreateGlobalVector(da,&global);CHKERRQ(ierr);
26   ierr = DMDAVecGetArrayDOF(da,global,&vglobal);CHKERRQ(ierr);
27   for (j=ystart; j<ystart+n; j++) {
28     for (i=xstart; i<xstart+m; i++) {
29       for (c=0; c<dof; c++) {
30         vglobal[j][i][c] = 100*j + 10*(i+1) + c;
31       }
32     }
33   }
34   ierr = DMDAVecRestoreArrayDOF(da,global,&vglobal);CHKERRQ(ierr);
35 
36   ierr = DMCreateLocalVector(da,&local);CHKERRQ(ierr);
37   ierr = DMGlobalToLocalBegin(da,global,INSERT_VALUES,local);CHKERRQ(ierr);
38   ierr = DMGlobalToLocalEnd(da,global,INSERT_VALUES,local);CHKERRQ(ierr);
39 
40   ierr = PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sview);CHKERRQ(ierr);
41   ierr = VecView(local,sview);CHKERRQ(ierr);
42   ierr = PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sview);CHKERRQ(ierr);
43   ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
44   ierr = VecView(global,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
45 
46   ierr = DMDestroy(&da);CHKERRQ(ierr);
47   ierr = VecDestroy(&local);CHKERRQ(ierr);
48   ierr = VecDestroy(&global);CHKERRQ(ierr);
49 
50   ierr = PetscFinalize();
51   return ierr;
52 }
53 
54 
55 
56 
57 
58 
59 
60 /*TEST
61 
62    test:
63 
64    test:
65       suffix: 2
66       nsize: 4
67       filter: grep -v "Vec Object"
68 
69 TEST*/
70