Actual source code: gr1.c

  1: #define PETSCDM_DLL

  3: /* 
  4:    Plots vectors obtained with DACreate1d()
  5: */

 7:  #include petscda.h

 11: /*@
 12:     DASetUniformCoordinates - Sets a DA coordinates to be a uniform grid

 14:   Collective on DA

 16:   Input Parameters:
 17: +  da - the distributed array object
 18: .  xmin,xmax - extremes in the x direction
 19: .  ymin,ymax - extremes in the y direction (use PETSC_NULL for 1 dimensional problems)
 20: -  zmin,zmax - extremes in the z direction (use PETSC_NULL for 1 or 2 dimensional problems)

 22:   Level: beginner

 24: .seealso: DASetCoordinates(), DAGetCoordinates(), DACreate1d(), DACreate2d(), DACreate3d()

 26: @*/
 27: PetscErrorCode  DASetUniformCoordinates(DA da,PetscReal xmin,PetscReal xmax,PetscReal ymin,PetscReal ymax,PetscReal zmin,PetscReal zmax)
 28: {
 29:   MPI_Comm       comm;
 30:   DA             cda;
 31:   DAPeriodicType periodic;
 32:   Vec            xcoor;
 33:   PetscScalar   *coors;
 34:   PetscReal      hx,hy,hz_;
 35:   PetscInt       i,j,k,M,N,P,istart,isize,jstart,jsize,kstart,ksize,dim,cnt;

 39:   if (xmax <= xmin) SETERRQ2(PETSC_ERR_ARG_INCOMP,"xmax must be larger than xmin %G %G",xmin,xmax);

 41:   PetscObjectGetComm((PetscObject)da,&comm);
 42:   DAGetInfo(da,&dim,&M,&N,&P,0,0,0,0,0,&periodic,0);
 43:   DAGetCorners(da,&istart,&jstart,&kstart,&isize,&jsize,&ksize);
 44:   DAGetCoordinateDA(da, &cda);
 45:   DACreateGlobalVector(cda, &xcoor);
 46:   if (dim == 1) {
 47:     if (periodic == DA_NONPERIODIC) hx = (xmax-xmin)/(M-1);
 48:     else                            hx = (xmax-xmin)/M;
 49:     VecGetArray(xcoor,&coors);
 50:     for (i=0; i<isize; i++) {
 51:       coors[i] = xmin + hx*(i+istart);
 52:     }
 53:     VecRestoreArray(xcoor,&coors);
 54:   } else if (dim == 2) {
 55:     if (ymax <= ymin) SETERRQ2(PETSC_ERR_ARG_INCOMP,"ymax must be larger than ymin %G %G",ymin,ymax);
 56:     if (DAXPeriodic(periodic)) hx = (xmax-xmin)/(M);
 57:     else                       hx = (xmax-xmin)/(M-1);
 58:     if (DAYPeriodic(periodic)) hy = (ymax-ymin)/(N);
 59:     else                       hy = (ymax-ymin)/(N-1);
 60:     VecGetArray(xcoor,&coors);
 61:     cnt  = 0;
 62:     for (j=0; j<jsize; j++) {
 63:       for (i=0; i<isize; i++) {
 64:         coors[cnt++] = xmin + hx*(i+istart);
 65:         coors[cnt++] = ymin + hy*(j+jstart);
 66:       }
 67:     }
 68:     VecRestoreArray(xcoor,&coors);
 69:   } else if (dim == 3) {
 70:     if (ymax <= ymin) SETERRQ2(PETSC_ERR_ARG_INCOMP,"ymax must be larger than ymin %G %G",ymin,ymax);
 71:     if (zmax <= zmin) SETERRQ2(PETSC_ERR_ARG_INCOMP,"zmax must be larger than zmin %G %G",zmin,zmax);
 72:     if (DAXPeriodic(periodic)) hx = (xmax-xmin)/(M);
 73:     else                       hx = (xmax-xmin)/(M-1);
 74:     if (DAYPeriodic(periodic)) hy = (ymax-ymin)/(N);
 75:     else                       hy = (ymax-ymin)/(N-1);
 76:     if (DAZPeriodic(periodic)) hz_ = (zmax-zmin)/(P);
 77:     else                       hz_ = (zmax-zmin)/(P-1);
 78:     VecGetArray(xcoor,&coors);
 79:     cnt  = 0;
 80:     for (k=0; k<ksize; k++) {
 81:       for (j=0; j<jsize; j++) {
 82:         for (i=0; i<isize; i++) {
 83:           coors[cnt++] = xmin + hx*(i+istart);
 84:           coors[cnt++] = ymin + hy*(j+jstart);
 85:           coors[cnt++] = zmin + hz_*(k+kstart);
 86:         }
 87:       }
 88:     }
 89:     VecRestoreArray(xcoor,&coors);
 90:   } else {
 91:     SETERRQ1(PETSC_ERR_SUP,"Cannot create uniform coordinates for this dimension %D\n",dim);
 92:   }
 93:   DASetCoordinates(da,xcoor);
 94:   PetscLogObjectParent(da,xcoor);
 95:   DADestroy(cda);
 96:   return(0);
 97: }

101: PetscErrorCode VecView_MPI_Draw_DA1d(Vec xin,PetscViewer v)
102: {
103:   DA             da;
105:   PetscMPIInt    rank,size,tag1,tag2;
106:   PetscInt       i,n,N,step,istart,isize,j;
107:   MPI_Status     status;
108:   PetscReal      coors[4],ymin,ymax,min,max,xmin,xmax,tmp,xgtmp;
109:   PetscScalar    *array,*xg;
110:   PetscDraw      draw;
111:   PetscTruth     isnull,showpoints;
112:   MPI_Comm       comm;
113:   PetscDrawAxis  axis;
114:   Vec            xcoor;
115:   DAPeriodicType periodic;

118:   PetscViewerDrawGetDraw(v,0,&draw);
119:   PetscDrawIsNull(draw,&isnull); if (isnull) return(0);

121:   PetscObjectQuery((PetscObject)xin,"DA",(PetscObject*)&da);
122:   if (!da) SETERRQ(PETSC_ERR_ARG_WRONG,"Vector not generated from a DA");

124:   PetscOptionsHasName(PETSC_NULL,"-draw_vec_mark_points",&showpoints);

126:   DAGetInfo(da,0,&N,0,0,0,0,0,&step,0,&periodic,0);
127:   DAGetCorners(da,&istart,0,0,&isize,0,0);
128:   VecGetArray(xin,&array);
129:   VecGetLocalSize(xin,&n);
130:   n    = n/step;

132:   /* get coordinates of nodes */
133:   DAGetCoordinates(da,&xcoor);
134:   if (!xcoor) {
135:     DASetUniformCoordinates(da,0.0,1.0,0.0,0.0,0.0,0.0);
136:     DAGetCoordinates(da,&xcoor);
137:   }
138:   VecGetArray(xcoor,&xg);

140:   PetscObjectGetComm((PetscObject)xin,&comm);
141:   MPI_Comm_size(comm,&size);
142:   MPI_Comm_rank(comm,&rank);

144:   /*
145:       Determine the min and max x coordinate in plot 
146:   */
147:   if (!rank) {
148:     xmin = PetscRealPart(xg[0]);
149:   }
150:   if (rank == size-1) {
151:     xmax = PetscRealPart(xg[n-1]);
152:   }
153:   MPI_Bcast(&xmin,1,MPIU_REAL,0,comm);
154:   MPI_Bcast(&xmax,1,MPIU_REAL,size-1,comm);

156:   for (j=0; j<step; j++) {
157:     PetscViewerDrawGetDraw(v,j,&draw);
158:     PetscDrawCheckResizedWindow(draw);

160:     /*
161:         Determine the min and max y coordinate in plot 
162:     */
163:     min = 1.e20; max = -1.e20;
164:     for (i=0; i<n; i++) {
165: #if defined(PETSC_USE_COMPLEX)
166:       if (PetscRealPart(array[j+i*step]) < min) min = PetscRealPart(array[j+i*step]);
167:       if (PetscRealPart(array[j+i*step]) > max) max = PetscRealPart(array[j+i*step]);
168: #else
169:       if (array[j+i*step] < min) min = array[j+i*step];
170:       if (array[j+i*step] > max) max = array[j+i*step];
171: #endif
172:     }
173:     if (min + 1.e-10 > max) {
174:       min -= 1.e-5;
175:       max += 1.e-5;
176:     }
177:     MPI_Reduce(&min,&ymin,1,MPIU_REAL,MPI_MIN,0,comm);
178:     MPI_Reduce(&max,&ymax,1,MPIU_REAL,MPI_MAX,0,comm);

180:     PetscDrawSynchronizedClear(draw);
181:     PetscViewerDrawGetDrawAxis(v,j,&axis);
182:     PetscLogObjectParent(draw,axis);
183:     if (!rank) {
184:       char *title;

186:       PetscDrawAxisSetLimits(axis,xmin,xmax,ymin,ymax);
187:       PetscDrawAxisDraw(axis);
188:       PetscDrawGetCoordinates(draw,coors,coors+1,coors+2,coors+3);
189:       DAGetFieldName(da,j,&title);
190:       if (title) {PetscDrawSetTitle(draw,title);}
191:     }
192:     MPI_Bcast(coors,4,MPIU_REAL,0,comm);
193:     if (rank) {
194:       PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);
195:     }

197:     /* draw local part of vector */
198:     PetscObjectGetNewTag((PetscObject)xin,&tag1);
199:     PetscObjectGetNewTag((PetscObject)xin,&tag2);
200:     if (rank < size-1) { /*send value to right */
201:       MPI_Send(&array[j+(n-1)*step],1,MPIU_REAL,rank+1,tag1,comm);
202:       MPI_Send(&xg[n-1],1,MPIU_REAL,rank+1,tag1,comm);
203:     }
204:     if (!rank && periodic && size > 1) { /* first processor sends first value to last */
205:       MPI_Send(&array[j],1,MPIU_REAL,size-1,tag2,comm);
206:     }

208:     for (i=1; i<n; i++) {
209: #if !defined(PETSC_USE_COMPLEX)
210:       PetscDrawLine(draw,xg[i-1],array[j+step*(i-1)],xg[i],array[j+step*i],
211:                       PETSC_DRAW_RED);
212: #else
213:       PetscDrawLine(draw,PetscRealPart(xg[i-1]),PetscRealPart(array[j+step*(i-1)]),
214:                       PetscRealPart(xg[i]),PetscRealPart(array[j+step*i]),PETSC_DRAW_RED);
215: #endif
216:       if (showpoints) {
217:         PetscDrawPoint(draw,PetscRealPart(xg[i-1]),PetscRealPart(array[j+step*(i-1)]),PETSC_DRAW_BLACK);
218:       }
219:     }
220:     if (rank) { /* receive value from left */
221:       MPI_Recv(&tmp,1,MPIU_REAL,rank-1,tag1,comm,&status);
222:       MPI_Recv(&xgtmp,1,MPIU_REAL,rank-1,tag1,comm,&status);
223: #if !defined(PETSC_USE_COMPLEX)
224:       PetscDrawLine(draw,xgtmp,tmp,xg[0],array[j],PETSC_DRAW_RED);
225: #else
226:       PetscDrawLine(draw,xgtmp,tmp,PetscRealPart(xg[0]),PetscRealPart(array[j]),
227:                       PETSC_DRAW_RED);
228: #endif
229:       if (showpoints) {
230:         PetscDrawPoint(draw,xgtmp,tmp,PETSC_DRAW_BLACK);
231:       }
232:     }
233:     if (rank == size-1 && periodic && size > 1) {
234:       MPI_Recv(&tmp,1,MPIU_REAL,0,tag2,comm,&status);
235: #if !defined(PETSC_USE_COMPLEX)
236:       PetscDrawLine(draw,xg[n-2],array[j+step*(n-1)],xg[n-1],tmp,PETSC_DRAW_RED);
237: #else
238:       PetscDrawLine(draw,PetscRealPart(xg[n-2]),PetscRealPart(array[j+step*(n-1)]),
239:                       PetscRealPart(xg[n-1]),tmp,PETSC_DRAW_RED);
240: #endif
241:       if (showpoints) {
242:         PetscDrawPoint(draw,PetscRealPart(xg[n-2]),PetscRealPart(array[j+step*(n-1)]),PETSC_DRAW_BLACK);
243:       }
244:     }
245:     PetscDrawSynchronizedFlush(draw);
246:     PetscDrawPause(draw);
247:   }
248:   VecRestoreArray(xcoor,&xg);
249:   VecRestoreArray(xin,&array);
250:   VecDestroy(xcoor);
251:   return(0);
252: }