Actual source code: ex5.c
petsc-3.8.4 2018-03-24
2: static char help[] = "Tests DMDAGetElements() and VecView() 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 M = 10,N = 8,ne,nc,i;
11: const PetscInt *e;
12: PetscErrorCode ierr;
13: PetscBool flg = PETSC_FALSE;
14: DM da;
15: PetscViewer viewer;
16: Vec local,global;
17: PetscScalar value;
18: DMBoundaryType bx = DM_BOUNDARY_NONE,by = DM_BOUNDARY_NONE;
19: DMDAStencilType stype = DMDA_STENCIL_BOX;
20: PetscScalar *lv;
22: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
23: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",300,0,300,300,&viewer);
25: /* Read options */
26: PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);
27: PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);
28: PetscOptionsGetBool(NULL,NULL,"-star_stencil",&flg,NULL);
29: if (flg) stype = DMDA_STENCIL_STAR;
31: /* Create distributed array and get vectors */
32: DMDACreate2d(PETSC_COMM_WORLD,bx,by,stype,M,N,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);
33: DMSetFromOptions(da);
34: DMSetUp(da);
35: DMCreateGlobalVector(da,&global);
36: DMCreateLocalVector(da,&local);
38: value = -3.0;
39: VecSet(global,value);
40: DMGlobalToLocalBegin(da,global,INSERT_VALUES,local);
41: DMGlobalToLocalEnd(da,global,INSERT_VALUES,local);
43: DMDASetElementType(da,DMDA_ELEMENT_P1);
44: DMDAGetElements(da,&ne,&nc,&e);
45: VecGetArray(local,&lv);
46: for (i=0; i<ne; i++) {
47: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"i %D e[3*i] %D %D %D\n",i,e[3*i],e[3*i+1],e[3*i+2]);
48: lv[e[3*i]] = i;
49: }
50: PetscSynchronizedFlush(PETSC_COMM_WORLD,stdout);
51: VecRestoreArray(local,&lv);
52: DMDARestoreElements(da,&ne,&nc,&e);
54: DMLocalToGlobalBegin(da,local,ADD_VALUES,global);
55: DMLocalToGlobalEnd(da,local,ADD_VALUES,global);
57: DMView(da,viewer);
58: VecView(global,viewer);
60: /* Free memory */
61: PetscViewerDestroy(&viewer);
62: VecDestroy(&local);
63: VecDestroy(&global);
64: DMDestroy(&da);
65: PetscFinalize();
66: return ierr;
67: }