Actual source code: ex40.c
petsc-3.8.4 2018-03-24
2: static char help[] = "Tests mirror boundary conditions in 2-d.\n\n";
4: #include <petscdm.h>
5: #include <petscdmda.h>
7: int main(int argc,char **argv)
8: {
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;
15: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
16: PetscOptionsGetInt(NULL,0,"-stencil_width",&stencil_width,0);
17: PetscOptionsGetInt(NULL,0,"-dof",&dof,0);
19: 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: DMSetFromOptions(da);
21: DMSetUp(da);
22: DMDAGetCorners(da,&xstart,&ystart,0,&m,&n,0);
24: DMCreateGlobalVector(da,&global);
25: 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: DMDAVecRestoreArrayDOF(da,global,&vglobal);
35: DMCreateLocalVector(da,&local);
36: DMGlobalToLocalBegin(da,global,INSERT_VALUES,local);
37: DMGlobalToLocalEnd(da,global,INSERT_VALUES,local);
39: PetscSequentialPhaseBegin(PETSC_COMM_WORLD,1);
40: VecView(local,PETSC_VIEWER_STDOUT_SELF);
41: PetscSequentialPhaseEnd(PETSC_COMM_WORLD,1);
42: VecView(global,PETSC_VIEWER_STDOUT_WORLD);
44: DMDestroy(&da);
45: VecDestroy(&local);
46: VecDestroy(&global);
48: PetscFinalize();
49: return ierr;
50: }