xref: /petsc/src/dm/tests/ex40.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
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   PetscFunctionBeginUser;
16   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
17   PetscCall(PetscOptionsGetInt(NULL,0,"-stencil_width",&stencil_width,0));
18   PetscCall(PetscOptionsGetInt(NULL,0,"-dof",&dof,0));
19 
20   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));
21   PetscCall(DMSetFromOptions(da));
22   PetscCall(DMSetUp(da));
23   PetscCall(DMDAGetCorners(da,&xstart,&ystart,0,&m,&n,0));
24 
25   PetscCall(DMCreateGlobalVector(da,&global));
26   PetscCall(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   PetscCall(DMDAVecRestoreArrayDOF(da,global,&vglobal));
35 
36   PetscCall(DMCreateLocalVector(da,&local));
37   PetscCall(DMGlobalToLocalBegin(da,global,INSERT_VALUES,local));
38   PetscCall(DMGlobalToLocalEnd(da,global,INSERT_VALUES,local));
39 
40   PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sview));
41   PetscCall(VecView(local,sview));
42   PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sview));
43   PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
44   PetscCall(VecView(global,PETSC_VIEWER_STDOUT_WORLD));
45 
46   PetscCall(DMDestroy(&da));
47   PetscCall(VecDestroy(&local));
48   PetscCall(VecDestroy(&global));
49 
50   PetscCall(PetscFinalize());
51   return 0;
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