Actual source code: ex41.c
petsc-3.11.4 2019-09-28
2: static char help[] = "Tests mirror boundary conditions in 3-d.\n\n";
4: #include <petscdm.h>
5: #include <petscdmda.h>
7: int main(int argc,char **argv)
8: {
10: PetscInt M = 2, N = 3, P = 4,stencil_width = 1, dof = 1,m,n,p,xstart,ystart,zstart,i,j,k,c;
11: DM da;
12: Vec global,local;
13: PetscScalar ****vglobal;
14: PetscViewer sview;
15: PetscScalar sum;
17: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
18: PetscOptionsGetInt(NULL,0,"-stencil_width",&stencil_width,0);
19: PetscOptionsGetInt(NULL,0,"-dof",&dof,0);
21: 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);
22: DMSetFromOptions(da);
23: DMSetUp(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 + c;
33: }
34: }
35: }
36: }
37: DMDAVecRestoreArrayDOF(da,global,&vglobal);
39: DMCreateLocalVector(da,&local);
40: DMGlobalToLocalBegin(da,global,ADD_VALUES,local);
41: DMGlobalToLocalEnd(da,global,ADD_VALUES,local);
43: VecSum(local,&sum);
44: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"sum %g\n",(double)PetscRealPart(sum));
45: PetscSynchronizedFlush(PETSC_COMM_WORLD,stdout);
46: PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sview);
47: VecView(local,sview);
48: PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sview);
49: PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);
51: DMDestroy(&da);
52: VecDestroy(&local);
53: VecDestroy(&global);
55: PetscFinalize();
56: return ierr;
57: }
59: /*TEST
61: test:
63: test:
64: suffix: 2
65: nsize: 3
67: TEST*/