Actual source code: ex41.c
petsc-3.6.1 2015-08-06
2: static char help[] = "Tests mirror boundary conditions in 3-d.\n\n";
4: #include <petscdm.h>
5: #include <petscdmda.h>
9: int main(int argc,char **argv)
10: {
12: PetscInt M = -2, N = -3, P = 4,stencil_width = 1, dof = 1,m,n,p,xstart,ystart,zstart,i,j,k,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: DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_MIRROR,DM_BOUNDARY_MIRROR,DM_BOUNDARY_MIRROR,DMDA_STENCIL_STAR,M,N,P,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,dof,stencil_width,NULL,NULL,NULL,&da);
24: DMDAGetCorners(da,&xstart,&ystart,&zstart,&m,&n,&p);
26: DMCreateGlobalVector(da,&global);
27: DMDAVecGetArrayDOF(da,global,&vglobal);
28: for (k=zstart; k<zstart+p; k++) {
29: for (j=ystart; j<ystart+n; j++) {
30: for (i=xstart; i<xstart+m; i++) {
31: for (c=0; c<dof; c++) {
32: vglobal[k][j][i][c] = 1000*k + 100*j + 10*(i+1) + c;
33: }
34: }
35: }
36: }
37: DMDAVecRestoreArrayDOF(da,global,&vglobal);
39: DMCreateLocalVector(da,&local);
40: DMGlobalToLocalBegin(da,global,INSERT_VALUES,local);
41: DMGlobalToLocalEnd(da,global,INSERT_VALUES,local);
43: PetscSequentialPhaseBegin(PETSC_COMM_WORLD,1);
44: VecView(local,PETSC_VIEWER_STDOUT_SELF);
45: PetscSequentialPhaseEnd(PETSC_COMM_WORLD,1);
46: VecView(global,PETSC_VIEWER_STDOUT_WORLD);
48: DMDestroy(&da);
49: VecDestroy(&local);
50: VecDestroy(&global);
52: PetscFinalize();
53: return 0;
54: }