xref: /petsc/src/dm/tests/ex40.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
1 
2 static char help[] = "Tests mirror boundary conditions in 2-d.\n\n";
3 
4 #include <petscdm.h>
5 #include <petscdmda.h>
6 
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   CHKERRQ(PetscOptionsGetInt(NULL,0,"-stencil_width",&stencil_width,0));
18   CHKERRQ(PetscOptionsGetInt(NULL,0,"-dof",&dof,0));
19 
20   CHKERRQ(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));
21   CHKERRQ(DMSetFromOptions(da));
22   CHKERRQ(DMSetUp(da));
23   CHKERRQ(DMDAGetCorners(da,&xstart,&ystart,0,&m,&n,0));
24 
25   CHKERRQ(DMCreateGlobalVector(da,&global));
26   CHKERRQ(DMDAVecGetArrayDOF(da,global,&vglobal));
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   CHKERRQ(DMDAVecRestoreArrayDOF(da,global,&vglobal));
35 
36   CHKERRQ(DMCreateLocalVector(da,&local));
37   CHKERRQ(DMGlobalToLocalBegin(da,global,INSERT_VALUES,local));
38   CHKERRQ(DMGlobalToLocalEnd(da,global,INSERT_VALUES,local));
39 
40   CHKERRQ(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sview));
41   CHKERRQ(VecView(local,sview));
42   CHKERRQ(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sview));
43   CHKERRQ(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
44   CHKERRQ(VecView(global,PETSC_VIEWER_STDOUT_WORLD));
45 
46   CHKERRQ(DMDestroy(&da));
47   CHKERRQ(VecDestroy(&local));
48   CHKERRQ(VecDestroy(&global));
49 
50   ierr = PetscFinalize();
51   return ierr;
52 }
53 
54 /*TEST
55 
56    test:
57 
58    test:
59       suffix: 2
60       nsize: 4
61       filter: grep -v "Vec Object"
62 
63 TEST*/
64