xref: /petsc/src/dm/tests/ex40.c (revision ffa8c5705e8ab2cf85ee1d14dbe507a6e2eb5283)
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   PetscInt       M = 8, N = 8,stencil_width = 1, dof = 1,m,n,xstart,ystart,i,j,c;
10   DM             da;
11   Vec            global,local;
12   PetscScalar    ***vglobal;
13   PetscViewer    sview;
14 
15   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
16   PetscCall(PetscOptionsGetInt(NULL,0,"-stencil_width",&stencil_width,0));
17   PetscCall(PetscOptionsGetInt(NULL,0,"-dof",&dof,0));
18 
19   PetscCall(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));
20   PetscCall(DMSetFromOptions(da));
21   PetscCall(DMSetUp(da));
22   PetscCall(DMDAGetCorners(da,&xstart,&ystart,0,&m,&n,0));
23 
24   PetscCall(DMCreateGlobalVector(da,&global));
25   PetscCall(DMDAVecGetArrayDOF(da,global,&vglobal));
26   for (j=ystart; j<ystart+n; j++) {
27     for (i=xstart; i<xstart+m; i++) {
28       for (c=0; c<dof; c++) {
29         vglobal[j][i][c] = 100*j + 10*(i+1) + c;
30       }
31     }
32   }
33   PetscCall(DMDAVecRestoreArrayDOF(da,global,&vglobal));
34 
35   PetscCall(DMCreateLocalVector(da,&local));
36   PetscCall(DMGlobalToLocalBegin(da,global,INSERT_VALUES,local));
37   PetscCall(DMGlobalToLocalEnd(da,global,INSERT_VALUES,local));
38 
39   PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sview));
40   PetscCall(VecView(local,sview));
41   PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sview));
42   PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
43   PetscCall(VecView(global,PETSC_VIEWER_STDOUT_WORLD));
44 
45   PetscCall(DMDestroy(&da));
46   PetscCall(VecDestroy(&local));
47   PetscCall(VecDestroy(&global));
48 
49   PetscCall(PetscFinalize());
50   return 0;
51 }
52 
53 /*TEST
54 
55    test:
56 
57    test:
58       suffix: 2
59       nsize: 4
60       filter: grep -v "Vec Object"
61 
62 TEST*/
63