Actual source code: ex39.c
petsc-3.7.3 2016-08-01
2: static char help[] = "Tests mirror boundary conditions in 1-d.\n\n";
4: #include <petscdm.h>
5: #include <petscdmda.h>
9: int main(int argc,char **argv)
10: {
12: PetscInt M = 6,stencil_width = 1, dof = 1,m,xstart,i,j;
13: DM da;
14: Vec global,local;
15: PetscScalar **vglobal;
17: PetscInitialize(&argc,&argv,(char*)0,help);
20: PetscOptionsGetInt(NULL,0,"-stencil_width",&stencil_width,0);
21: PetscOptionsGetInt(NULL,0,"-dof",&dof,0);
23: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_MIRROR,M,dof,stencil_width,NULL,&da);
24: DMDAGetCorners(da,&xstart,0,0,&m,0,0);
26: DMCreateGlobalVector(da,&global);
27: DMDAVecGetArrayDOF(da,global,&vglobal);
28: for (i=xstart; i<xstart+m; i++) {
29: for (j=0; j<dof; j++) {
30: vglobal[i][j] = 100*(i+1) + j;
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 0;
50: }