Actual source code: ex40.c
petsc-3.6.1 2015-08-06
2: static char help[] = "Tests mirror boundary conditions in 2-d.\n\n";
4: #include <petscdm.h>
5: #include <petscdmda.h>
9: int main(int argc,char **argv)
10: {
12: PetscInt M = -8, N = -8,stencil_width = 1, dof = 1,m,n,xstart,ystart,i,j,c;
13: DM da;
14: Vec global,local;
15: PetscScalar ***vglobal;
17: PetscInitialize(&argc,&argv,(char*)0,help);
20: PetscOptionsGetInt(0,"-stencil_width",&stencil_width,0);
21: PetscOptionsGetInt(0,"-dof",&dof,0);
23: 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);
24: DMDAGetCorners(da,&xstart,&ystart,0,&m,&n,0);
26: DMCreateGlobalVector(da,&global);
27: DMDAVecGetArrayDOF(da,global,&vglobal);
28: for (j=ystart; j<ystart+n; j++) {
29: for (i=xstart; i<xstart+m; i++) {
30: for (c=0; c<dof; c++) {
31: vglobal[j][i][c] = 100*j + 10*(i+1) + c;
32: }
33: }
34: }
35: DMDAVecRestoreArrayDOF(da,global,&vglobal);
37: DMCreateLocalVector(da,&local);
38: DMGlobalToLocalBegin(da,global,INSERT_VALUES,local);
39: DMGlobalToLocalEnd(da,global,INSERT_VALUES,local);
41: PetscSequentialPhaseBegin(PETSC_COMM_WORLD,1);
42: VecView(local,PETSC_VIEWER_STDOUT_SELF);
43: PetscSequentialPhaseEnd(PETSC_COMM_WORLD,1);
44: VecView(global,PETSC_VIEWER_STDOUT_WORLD);
46: DMDestroy(&da);
47: VecDestroy(&local);
48: VecDestroy(&global);
50: PetscFinalize();
51: return 0;
52: }