Actual source code: ex4.c
1:
2: static char help[] = "Tests various 2-dimensional DA routines.\n\n";
4: #include petscda.h
8: int main(int argc,char **argv)
9: {
10: PetscMPIInt rank;
12: PetscInt M = 10,N = 8,m = PETSC_DECIDE;
13: PetscInt s=2,w=2,n = PETSC_DECIDE,nloc,l,i,j,kk;
14: PetscInt Xs,Xm,Ys,Ym,iloc,*iglobal,*ltog;
15: PetscInt *lx = PETSC_NULL,*ly = PETSC_NULL;
16: PetscTruth testorder = PETSC_FALSE,flg;
17: DAPeriodicType wrap = DA_NONPERIODIC;
18: DA da;
19: PetscViewer viewer;
20: Vec local,global;
21: PetscScalar value;
22: DAStencilType st = DA_STENCIL_BOX;
23: AO ao;
24:
25: PetscInitialize(&argc,&argv,(char*)0,help);
26: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",300,0,400,400,&viewer);
27:
28: /* Readoptions */
29: PetscOptionsGetInt(PETSC_NULL,"-M",&M,PETSC_NULL);
30: PetscOptionsGetInt(PETSC_NULL,"-N",&N,PETSC_NULL);
31: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
32: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
33: PetscOptionsGetInt(PETSC_NULL,"-s",&s,PETSC_NULL);
34: PetscOptionsGetInt(PETSC_NULL,"-w",&w,PETSC_NULL);
35: flg = PETSC_FALSE;
36: PetscOptionsGetTruth(PETSC_NULL,"-xwrap",&flg,PETSC_NULL); if (flg) wrap = DA_XPERIODIC;
37: flg = PETSC_FALSE;
38: PetscOptionsGetTruth(PETSC_NULL,"-ywrap",&flg,PETSC_NULL); if (flg) wrap = DA_YPERIODIC;
39: flg = PETSC_FALSE;
40: PetscOptionsGetTruth(PETSC_NULL,"-xywrap",&flg,PETSC_NULL); if (flg) wrap = DA_XYPERIODIC;
41: flg = PETSC_FALSE;
42: PetscOptionsGetTruth(PETSC_NULL,"-star",&flg,PETSC_NULL); if (flg) st = DA_STENCIL_STAR;
43: flg = PETSC_FALSE;
44: PetscOptionsGetTruth(PETSC_NULL,"-testorder",&testorder,PETSC_NULL);
45: /*
46: Test putting two nodes in x and y on each processor, exact last processor
47: in x and y gets the rest.
48: */
49: flg = PETSC_FALSE;
50: PetscOptionsGetTruth(PETSC_NULL,"-distribute",&flg,PETSC_NULL);
51: if (flg) {
52: if (m == PETSC_DECIDE) SETERRQ(1,"Must set -m option with -distribute option");
53: PetscMalloc(m*sizeof(PetscInt),&lx);
54: for (i=0; i<m-1; i++) { lx[i] = 4;}
55: lx[m-1] = M - 4*(m-1);
56: if (n == PETSC_DECIDE) SETERRQ(1,"Must set -n option with -distribute option");
57: PetscMalloc(n*sizeof(PetscInt),&ly);
58: for (i=0; i<n-1; i++) { ly[i] = 2;}
59: ly[n-1] = N - 2*(n-1);
60: }
63: /* Create distributed array and get vectors */
64: DACreate2d(PETSC_COMM_WORLD,wrap,st,M,N,m,n,w,s,lx,ly,&da);
65: if (lx) {
66: PetscFree(lx);
67: PetscFree(ly);
68: }
70: DAView(da,viewer);
71: DACreateGlobalVector(da,&global);
72: DACreateLocalVector(da,&local);
74: /* Set global vector; send ghost points to local vectors */
75: value = 1;
76: VecSet(global,value);
77: DAGlobalToLocalBegin(da,global,INSERT_VALUES,local);
78: DAGlobalToLocalEnd(da,global,INSERT_VALUES,local);
80: /* Scale local vectors according to processor rank; pass to global vector */
81: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
82: value = rank;
83: VecScale(local,value);
84: DALocalToGlobal(da,local,INSERT_VALUES,global);
86: if (!testorder) { /* turn off printing when testing ordering mappings */
87: PetscPrintf (PETSC_COMM_WORLD,"\nGlobal Vectors:\n");
88: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_NATIVE);
89: VecView(global,PETSC_VIEWER_STDOUT_WORLD);
90: PetscPrintf (PETSC_COMM_WORLD,"\n\n");
91: }
93: /* Send ghost points to local vectors */
94: DAGlobalToLocalBegin(da,global,INSERT_VALUES,local);
95: DAGlobalToLocalEnd(da,global,INSERT_VALUES,local);
97: flg = PETSC_FALSE;
98: PetscOptionsGetTruth(PETSC_NULL,"-local_print",&flg,PETSC_NULL);
99: if (flg) {
100: PetscViewer sviewer;
101: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"\nLocal Vector: processor %d\n",rank);
102: PetscViewerGetSingleton(PETSC_VIEWER_STDOUT_WORLD,&sviewer);
103: VecView(local,sviewer);
104: PetscViewerRestoreSingleton(PETSC_VIEWER_STDOUT_WORLD,&sviewer);
105: }
107: /* Tests mappings betweeen application/PETSc orderings */
108: if (testorder) {
109: DAGetGhostCorners(da,&Xs,&Ys,PETSC_NULL,&Xm,&Ym,PETSC_NULL);
110: DAGetGlobalIndices(da,&nloc,<og);
111: DAGetAO(da,&ao);
112: PetscMalloc(nloc*sizeof(PetscInt),&iglobal);
114: /* Set iglobal to be global indices for each processor's local and ghost nodes,
115: using the DA ordering of grid points */
116: kk = 0;
117: for (j=Ys; j<Ys+Ym; j++) {
118: for (i=Xs; i<Xs+Xm; i++) {
119: iloc = w*((j-Ys)*Xm + i-Xs);
120: for (l=0; l<w; l++) {
121: iglobal[kk++] = ltog[iloc+l];
122: }
123: }
124: }
126: /* Map this to the application ordering (which for DAs is just the natural ordering
127: that would be used for 1 processor, numbering most rapidly by x, then y) */
128: AOPetscToApplication(ao,nloc,iglobal);
130: /* Then map the application ordering back to the PETSc DA ordering */
131: AOApplicationToPetsc(ao,nloc,iglobal);
133: /* Verify the mappings */
134: kk=0;
135: for (j=Ys; j<Ys+Ym; j++) {
136: for (i=Xs; i<Xs+Xm; i++) {
137: iloc = w*((j-Ys)*Xm + i-Xs);
138: for (l=0; l<w; l++) {
139: if (iglobal[kk] != ltog[iloc+l]) {
140: PetscFPrintf(PETSC_COMM_SELF,stdout,"[%d] Problem with mapping: j=%D, i=%D, l=%D, petsc1=%D, petsc2=%D\n",
141: rank,j,i,l,ltog[iloc+l],iglobal[kk]);}
142: kk++;
143: }
144: }
145: }
146: PetscFree(iglobal);
147: }
149: /* Free memory */
150: PetscViewerDestroy(viewer);
151: VecDestroy(local);
152: VecDestroy(global);
153: DADestroy(da);
155: PetscFinalize();
156: return 0;
157: }