Actual source code: ex2.c
petsc-3.8.4 2018-03-24
2: static char help[] = "Tests DMDAGlobalToNaturalAllCreate() using contour plotting for 2d DMDAs.\n\n";
4: #include <petscdm.h>
5: #include <petscdmda.h>
6: #include <petscdraw.h>
8: int main(int argc,char **argv)
9: {
10: PetscInt i,j,M = 10,N = 8,m = PETSC_DECIDE,n = PETSC_DECIDE;
11: PetscMPIInt rank;
12: PetscErrorCode ierr;
13: PetscBool flg = PETSC_FALSE;
14: DM da;
15: PetscViewer viewer;
16: Vec localall,global;
17: PetscScalar value,*vlocal;
18: DMBoundaryType bx = DM_BOUNDARY_NONE,by = DM_BOUNDARY_NONE;
19: DMDAStencilType stype = DMDA_STENCIL_BOX;
20: VecScatter tolocalall,fromlocalall;
21: PetscInt start,end;
23: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
24: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",300,0,300,300,&viewer);
26: /* Read options */
27: PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);
28: PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);
29: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
30: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
31: PetscOptionsGetBool(NULL,NULL,"-star_stencil",&flg,NULL);
32: if (flg) stype = DMDA_STENCIL_STAR;
34: /* Create distributed array and get vectors */
35: DMDACreate2d(PETSC_COMM_WORLD,bx,by,stype,M,N,m,n,1,1,NULL,NULL,&da);
36: DMSetFromOptions(da);
37: DMSetUp(da);
38:
39: DMCreateGlobalVector(da,&global);
40: VecCreateSeq(PETSC_COMM_SELF,M*N,&localall);
42: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
43: VecGetOwnershipRange(global,&start,&end);
44: for (i=start; i<end; i++) {
45: value = 5.0*rank;
46: VecSetValues(global,1,&i,&value,INSERT_VALUES);
47: }
48: VecView(global,viewer);
50: /*
51: Create Scatter from global DMDA parallel vector to local vector that
52: contains all entries
53: */
54: DMDAGlobalToNaturalAllCreate(da,&tolocalall);
55: DMDANaturalAllToGlobalCreate(da,&fromlocalall);
57: VecScatterBegin(tolocalall,global,localall,INSERT_VALUES,SCATTER_FORWARD);
58: VecScatterEnd(tolocalall,global,localall,INSERT_VALUES,SCATTER_FORWARD);
60: VecGetArray(localall,&vlocal);
61: for (j=0; j<N; j++) {
62: for (i=0; i<M; i++) {
63: *vlocal++ += i + j*M;
64: }
65: }
66: VecRestoreArray(localall,&vlocal);
68: /* scatter back to global vector */
69: VecScatterBegin(fromlocalall,localall,global,INSERT_VALUES,SCATTER_FORWARD);
70: VecScatterEnd(fromlocalall,localall,global,INSERT_VALUES,SCATTER_FORWARD);
72: VecView(global,viewer);
74: /* Free memory */
75: VecScatterDestroy(&tolocalall);
76: VecScatterDestroy(&fromlocalall);
77: PetscViewerDestroy(&viewer);
78: VecDestroy(&localall);
79: VecDestroy(&global);
80: DMDestroy(&da);
81: PetscFinalize();
82: return ierr;
83: }