1 2 static char help[] = "Tests mirror boundary conditions in 3-d.\n\n"; 3 4 #include <petscdm.h> 5 #include <petscdmda.h> 6 7 int main(int argc,char **argv) 8 { 9 PetscInt M = 2, N = 3, P = 4,stencil_width = 1, dof = 1,m,n,p,xstart,ystart,zstart,i,j,k,c; 10 DM da; 11 Vec global,local; 12 PetscScalar ****vglobal; 13 PetscViewer sview; 14 PetscScalar sum; 15 16 PetscFunctionBeginUser; 17 PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 18 PetscCall(PetscOptionsGetInt(NULL,0,"-stencil_width",&stencil_width,0)); 19 PetscCall(PetscOptionsGetInt(NULL,0,"-dof",&dof,0)); 20 21 PetscCall(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 PetscCall(DMSetFromOptions(da)); 23 PetscCall(DMSetUp(da)); 24 PetscCall(DMDAGetCorners(da,&xstart,&ystart,&zstart,&m,&n,&p)); 25 26 PetscCall(DMCreateGlobalVector(da,&global)); 27 PetscCall(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 PetscCall(DMDAVecRestoreArrayDOF(da,global,&vglobal)); 38 39 PetscCall(DMCreateLocalVector(da,&local)); 40 PetscCall(DMGlobalToLocalBegin(da,global,ADD_VALUES,local)); 41 PetscCall(DMGlobalToLocalEnd(da,global,ADD_VALUES,local)); 42 43 PetscCall(VecSum(local,&sum)); 44 PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"sum %g\n",(double)PetscRealPart(sum))); 45 PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,stdout)); 46 PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sview)); 47 PetscCall(VecView(local,sview)); 48 PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sview)); 49 PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 50 51 PetscCall(DMDestroy(&da)); 52 PetscCall(VecDestroy(&local)); 53 PetscCall(VecDestroy(&global)); 54 55 PetscCall(PetscFinalize()); 56 return 0; 57 } 58 59 /*TEST 60 61 test: 62 63 test: 64 suffix: 2 65 nsize: 3 66 67 TEST*/ 68