Actual source code: gr2.c

  1: /* 
  2:    Plots vectors obtained with DACreate2d()
  3: */

 5:  #include src/dm/da/daimpl.h
 6:  #include vecimpl.h

  8: #if defined(PETSC_HAVE_PNETCDF)
 10: #include "pnetcdf.h"
 12: #endif


 15: /*
 16:         The data that is passed into the graphics callback
 17: */
 18: typedef struct {
 19:   PetscInt          m,n,step,k;
 20:   PetscReal    min,max,scale;
 21:   PetscScalar  *xy,*v;
 22:   PetscTruth   showgrid;
 23: } ZoomCtx;

 25: /*
 26:        This does the drawing for one particular field 
 27:     in one particular set of coordinates. It is a callback
 28:     called from PetscDrawZoom()
 29: */
 32: PetscErrorCode VecView_MPI_Draw_DA2d_Zoom(PetscDraw draw,void *ctx)
 33: {
 34:   ZoomCtx     *zctx = (ZoomCtx*)ctx;
 36:   PetscInt m,n,i,j,k,step,id,c1,c2,c3,c4;
 37:   PetscReal   s,min,x1,x2,x3,x4,y_1,y2,y3,y4;
 38:   PetscScalar *v,*xy;

 41:   m    = zctx->m;
 42:   n    = zctx->n;
 43:   step = zctx->step;
 44:   k    = zctx->k;
 45:   v    = zctx->v;
 46:   xy   = zctx->xy;
 47:   s    = zctx->scale;
 48:   min  = zctx->min;
 49: 
 50:   /* PetscDraw the contour plot patch */
 51:   for (j=0; j<n-1; j++) {
 52:     for (i=0; i<m-1; i++) {
 53: #if !defined(PETSC_USE_COMPLEX)
 54:       id = i+j*m;    x1 = xy[2*id];y_1 = xy[2*id+1];c1 = (int)(PETSC_DRAW_BASIC_COLORS+s*(v[k+step*id]-min));
 55:       id = i+j*m+1;  x2 = xy[2*id];y2  = y_1;       c2 = (int)(PETSC_DRAW_BASIC_COLORS+s*(v[k+step*id]-min));
 56:       id = i+j*m+1+m;x3 = x2;      y3  = xy[2*id+1];c3 = (int)(PETSC_DRAW_BASIC_COLORS+s*(v[k+step*id]-min));
 57:       id = i+j*m+m;  x4 = x1;      y4  = y3;        c4 = (int)(PETSC_DRAW_BASIC_COLORS+s*(v[k+step*id]-min));
 58: #else
 59:       id = i+j*m;    x1 = PetscRealPart(xy[2*id]);y_1 = PetscRealPart(xy[2*id+1]);c1 = (int)(PETSC_DRAW_BASIC_COLORS+s*(PetscRealPart(v[k+step*id])-min));
 60:       id = i+j*m+1;  x2 = PetscRealPart(xy[2*id]);y2  = y_1;       c2 = (int)(PETSC_DRAW_BASIC_COLORS+s*(PetscRealPart(v[k+step*id])-min));
 61:       id = i+j*m+1+m;x3 = x2;      y3  = PetscRealPart(xy[2*id+1]);c3 = (int)(PETSC_DRAW_BASIC_COLORS+s*(PetscRealPart(v[k+step*id])-min));
 62:       id = i+j*m+m;  x4 = x1;      y4  = y3;        c4 = (int)(PETSC_DRAW_BASIC_COLORS+s*(PetscRealPart(v[k+step*id])-min));
 63: #endif
 64:       PetscDrawTriangle(draw,x1,y_1,x2,y2,x3,y3,c1,c2,c3);
 65:       PetscDrawTriangle(draw,x1,y_1,x3,y3,x4,y4,c1,c3,c4);
 66:       if (zctx->showgrid) {
 67:         PetscDrawLine(draw,x1,y_1,x2,y2,PETSC_DRAW_BLACK);
 68:         PetscDrawLine(draw,x2,y2,x3,y3,PETSC_DRAW_BLACK);
 69:         PetscDrawLine(draw,x3,y3,x4,y4,PETSC_DRAW_BLACK);
 70:         PetscDrawLine(draw,x4,y4,x1,y_1,PETSC_DRAW_BLACK);
 71:       }
 72:     }
 73:   }
 74:   return(0);
 75: }

 79: PetscErrorCode VecView_MPI_Draw_DA2d(Vec xin,PetscViewer viewer)
 80: {
 81:   DA                 da,dac,dag;
 82:   PetscErrorCode     ierr;
 83:   PetscMPIInt        rank;
 84:   PetscInt           igstart,N,s,M,istart,isize,jgstart,*lx,*ly,w;
 85:   PetscReal          coors[4],ymin,ymax,xmin,xmax;
 86:   PetscDraw          draw,popup;
 87:   PetscTruth         isnull,useports;
 88:   MPI_Comm           comm;
 89:   Vec                xlocal,xcoor,xcoorl;
 90:   DAPeriodicType     periodic;
 91:   DAStencilType      st;
 92:   ZoomCtx            zctx;
 93:   PetscDrawViewPorts *ports;
 94:   PetscViewerFormat  format;

 97:   PetscViewerDrawGetDraw(viewer,0,&draw);
 98:   PetscDrawIsNull(draw,&isnull); if (isnull) return(0);

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

103:   PetscObjectGetComm((PetscObject)xin,&comm);
104:   MPI_Comm_rank(comm,&rank);

106:   DAGetInfo(da,0,&M,&N,0,&zctx.m,&zctx.n,0,&w,&s,&periodic,&st);
107:   DAGetOwnershipRange(da,&lx,&ly,PETSC_NULL);

109:   /* 
110:         Obtain a sequential vector that is going to contain the local values plus ONE layer of 
111:      ghosted values to draw the graphics from. We also need its corresponding DA (dac) that will
112:      update the local values pluse ONE layer of ghost values. 
113:   */
114:   PetscObjectQuery((PetscObject)da,"GraphicsGhosted",(PetscObject*)&xlocal);
115:   if (!xlocal) {
116:     if (periodic != DA_NONPERIODIC || s != 1 || st != DA_STENCIL_BOX) {
117:       /* 
118:          if original da is not of stencil width one, or periodic or not a box stencil then
119:          create a special DA to handle one level of ghost points for graphics
120:       */
121:       DACreate2d(comm,DA_NONPERIODIC,DA_STENCIL_BOX,M,N,zctx.m,zctx.n,w,1,lx,ly,&dac);
122:       PetscLogInfo(da,"VecView_MPI_Draw_DA2d:Creating auxilary DA for managing graphics ghost points\n");
123:     } else {
124:       /* otherwise we can use the da we already have */
125:       dac = da;
126:     }
127:     /* create local vector for holding ghosted values used in graphics */
128:     DACreateLocalVector(dac,&xlocal);
129:     if (dac != da) {
130:       /* don't keep any public reference of this DA, is is only available through xlocal */
131:       DADestroy(dac);
132:     } else {
133:       /* remove association between xlocal and da, because below we compose in the opposite
134:          direction and if we left this connect we'd get a loop, so the objects could 
135:          never be destroyed */
136:       PetscObjectCompose((PetscObject)xlocal,"DA",0);
137:     }
138:     PetscObjectCompose((PetscObject)da,"GraphicsGhosted",(PetscObject)xlocal);
139:     PetscObjectDereference((PetscObject)xlocal);
140:   } else {
141:     if (periodic == DA_NONPERIODIC && s == 1 && st == DA_STENCIL_BOX) {
142:       dac = da;
143:     } else {
144:       PetscObjectQuery((PetscObject)xlocal,"DA",(PetscObject*)&dac);
145:     }
146:   }

148:   /*
149:       Get local (ghosted) values of vector
150:   */
151:   DAGlobalToLocalBegin(dac,xin,INSERT_VALUES,xlocal);
152:   DAGlobalToLocalEnd(dac,xin,INSERT_VALUES,xlocal);
153:   VecGetArray(xlocal,&zctx.v);

155:   /* get coordinates of nodes */
156:   DAGetCoordinates(da,&xcoor);
157:   if (!xcoor) {
158:     DASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,0.0);
159:     DAGetCoordinates(da,&xcoor);
160:   }

162:   /*
163:       Determine the min and max  coordinates in plot 
164:   */
165:   VecStrideMin(xcoor,0,PETSC_NULL,&xmin);
166:   VecStrideMax(xcoor,0,PETSC_NULL,&xmax);
167:   VecStrideMin(xcoor,1,PETSC_NULL,&ymin);
168:   VecStrideMax(xcoor,1,PETSC_NULL,&ymax);
169:   coors[0] = xmin - .05*(xmax- xmin); coors[2] = xmax + .05*(xmax - xmin);
170:   coors[1] = ymin - .05*(ymax- ymin); coors[3] = ymax + .05*(ymax - ymin);
171:   PetscLogInfo(da,"VecView_MPI_Draw_DA2d:Preparing DA 2d contour plot coordinates %g %g %g %g\n",coors[0],coors[1],coors[2],coors[3]);

173:   /*
174:        get local ghosted version of coordinates 
175:   */
176:   PetscObjectQuery((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject*)&xcoorl);
177:   if (!xcoorl) {
178:     /* create DA to get local version of graphics */
179:     DACreate2d(comm,DA_NONPERIODIC,DA_STENCIL_BOX,M,N,zctx.m,zctx.n,2,1,lx,ly,&dag);
180:     PetscLogInfo(dag,"VecView_MPI_Draw_DA2d:Creating auxilary DA for managing graphics coordinates ghost points\n");
181:     DACreateLocalVector(dag,&xcoorl);
182:     PetscObjectCompose((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject)xcoorl);
183:     DADestroy(dag);/* dereference dag */
184:     PetscObjectDereference((PetscObject)xcoorl);
185:   } else {
186:     PetscObjectQuery((PetscObject)xcoorl,"DA",(PetscObject*)&dag);
187:   }
188:   DAGlobalToLocalBegin(dag,xcoor,INSERT_VALUES,xcoorl);
189:   DAGlobalToLocalEnd(dag,xcoor,INSERT_VALUES,xcoorl);
190:   VecGetArray(xcoorl,&zctx.xy);
191: 
192:   /*
193:         Get information about size of area each processor must do graphics for
194:   */
195:   DAGetInfo(dac,0,&M,&N,0,0,0,0,&zctx.step,0,&periodic,0);
196:   DAGetGhostCorners(dac,&igstart,&jgstart,0,&zctx.m,&zctx.n,0);
197:   DAGetCorners(dac,&istart,0,0,&isize,0,0);

199:   PetscOptionsHasName(PETSC_NULL,"-draw_contour_grid",&zctx.showgrid);

201:   PetscViewerGetFormat(viewer,&format);
202:   PetscOptionsHasName(PETSC_NULL,"-draw_ports",&useports);
203:   if (useports || format == PETSC_VIEWER_DRAW_PORTS){
204:     PetscDrawSynchronizedClear(draw);
205:     PetscDrawViewPortsCreate(draw,zctx.step,&ports);
206:   }
207:   /*
208:      Loop over each field; drawing each in a different window
209:   */
210:   for (zctx.k=0; zctx.k<zctx.step; zctx.k++) {
211:     if (useports) {
212:       PetscDrawViewPortsSet(ports,zctx.k);
213:     } else {
214:       PetscViewerDrawGetDraw(viewer,zctx.k,&draw);
215:       PetscDrawSynchronizedClear(draw);
216:     }

218:     /*
219:         Determine the min and max color in plot 
220:     */
221:     VecStrideMin(xin,zctx.k,PETSC_NULL,&zctx.min);
222:     VecStrideMax(xin,zctx.k,PETSC_NULL,&zctx.max);
223:     if (zctx.min == zctx.max) {
224:       zctx.min -= 1.e-12;
225:       zctx.max += 1.e-12;
226:     }

228:     if (!rank) {
229:       char *title;

231:       DAGetFieldName(da,zctx.k,&title);
232:       if (title) {
233:         PetscDrawSetTitle(draw,title);
234:       }
235:     }
236:     PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);
237:     PetscLogInfo(da,"VecView_MPI_Draw_DA2d:DA 2d contour plot min %g max %g\n",zctx.min,zctx.max);

239:     PetscDrawGetPopup(draw,&popup);
240:     if (popup) {PetscDrawScalePopup(popup,zctx.min,zctx.max);}

242:     zctx.scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/(zctx.max - zctx.min);

244:     PetscDrawZoom(draw,VecView_MPI_Draw_DA2d_Zoom,&zctx);
245:   }
246:   if (useports){
247:     PetscDrawViewPortsDestroy(ports);
248:   }

250:   VecRestoreArray(xcoorl,&zctx.xy);
251:   VecRestoreArray(xlocal,&zctx.v);
252:   return(0);
253: }

255: EXTERN PetscErrorCode VecView_MPI_HDF4_Ex(Vec X, PetscViewer viewer, PetscInt d, PetscInt *dims);

259: PetscErrorCode VecView_MPI_HDF4_DA2d(Vec xin,PetscViewer viewer)
260: {
261: #if defined(PETSC_HAVE_HDF4) && !defined(PETSC_USE_COMPLEX)
263:   PetscInt dims[2];
264:   DA  da;
265:   Vec natural;


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

272:   dims[0] = da->M;
273:   dims[1] = da->N;

275:   DACreateNaturalVector(da,&natural);
276:   DAGlobalToNaturalBegin(da,xin,INSERT_VALUES,natural);
277:   DAGlobalToNaturalEnd(da,xin,INSERT_VALUES,natural);
278:   VecView_MPI_HDF4_Ex(natural, viewer, 2, dims);
279:   VecDestroy(natural);

281:   return(0);
282: #else /* !defined(PETSC_HAVE_HDF4) */
284:   SETERRQ(PETSC_ERR_SUP_SYS,"Build PETSc with HDF4 to use this viewer");
285: #endif    
286: }

290: PetscErrorCode VecView_MPI_Netcdf_DA(Vec xin,PetscViewer viewer)
291: {
292: #if defined(PETSC_HAVE_PNETCDF)
294:   PetscInt ncid,xstart,xdim_num=1;
295:   PetscInt            i,j,len,dim,m,n,p,dof,swidth,M,N,P;
296:   PetscInt            xin_dim,xin_id,xin_n,xin_N,xyz_dim,xyz_id,xyz_n,xyz_N;
297:   PetscInt            *lx,*ly,*lz;
298:   PetscScalar    *xarray;
299:   DA             da,dac;
300:   Vec            natural,xyz;
301:   DAStencilType  stencil;
302:   DAPeriodicType periodic;
303:   MPI_Comm       comm;

306:   PetscObjectGetComm((PetscObject)xin,&comm);
307:   PetscObjectQuery((PetscObject)xin,"DA",(PetscObject*)&da);
308:   if (!da) SETERRQ(PETSC_ERR_ARG_WRONG,,"Vector not generated from a DA");
309:   DAGetInfo(da,&dim,&m,&n,&p,&M,&N,&P,&dof,&swidth,&periodic,&stencil);

311:   /* create the appropriate DA to map the coordinates to natural ordering */
312:   DAGetOwnershipRange(da,&lx,&ly,&lz);
313:   if (dim == 1) {
314:     DACreate1d(comm,DA_NONPERIODIC,m,dim,0,lx,&dac);
315:   } else if (dim == 2) {
316:     DACreate2d(comm,DA_NONPERIODIC,DA_STENCIL_BOX,m,n,M,N,dim,0,lx,ly,&dac);
317:   } else if (dim == 3) {
318:     DACreate3d(comm,DA_NONPERIODIC,DA_STENCIL_BOX,m,n,p,M,N,P,dim,0,lx,ly,lz,&dac);
319:   } else {
320:     SETERRQ1(PETSC_ERR_ARG_CORRUPT,"Dimension is not 1 2 or 3: %D\n",dim);
321:   }
322:   DACreateNaturalVector(dac,&xyz);
323:   PetscObjectSetOptionsPrefix((PetscObject)xyz,"coor_");
324:   DAGlobalToNaturalBegin(dac,da->coordinates,INSERT_VALUES,xyz);
325:   DAGlobalToNaturalEnd(dac,da->coordinates,INSERT_VALUES,xyz);
326:   /* Create the DA vector in natural ordering */
327:   DACreateNaturalVector(da,&natural);
328:   DAGlobalToNaturalBegin(da,xin,INSERT_VALUES,natural);
329:   DAGlobalToNaturalEnd(da,xin,INSERT_VALUES,natural);
330:   /* Write the netCDF dataset */
331:   PetscViewerNetcdfGetID(viewer,&ncid);
332:   if (ncid < 0) SETERRQ(PETSC_ERR_ORDER,"First call PetscViewerNetcdfOpen to create NetCDF dataset");
333:   /* define dimensions */
334:   VecGetSize(xin,&xin_N);
335:   VecGetLocalSize(xin,&xin_n);
336:   ncmpi_def_dim(ncid,"PETSc_DA_Vector_Global_Size",xin_N,&xin_dim);
337:   VecGetSize(xyz,&xyz_N);
338:   VecGetLocalSize(xyz,&xyz_n);
339:   ncmpi_def_dim(ncid,"PETSc_DA_Coordinate_Vector_Global_Size",xyz_N,&xyz_dim);
340:   /* define variables */
341:   ncmpi_def_var(ncid,"PETSc_DA_Vector",NC_DOUBLE,xdim_num,&xin_dim,&xin_id);
342:   ncmpi_def_var(ncid,"PETSc_DA_Coordinate_Vector",NC_DOUBLE,xdim_num,&xyz_dim,&xyz_id);
343:   /* leave define mode */
344:   ncmpi_enddef(ncid);
345:   /* store the vector */
346:   VecGetArray(xin,&xarray);
347:   VecGetOwnershipRange(xin,&xstart,PETSC_NULL);
348:   ncmpi_put_vara_double_all(ncid,xin_id,(const size_t*)&xstart,(const size_t*)&xin_n,xarray);
349:   VecRestoreArray(xin,&xarray);
350:   /* store the coordinate vector */
351:   VecGetArray(xyz,&xarray);
352:   VecGetOwnershipRange(xyz,&xstart,PETSC_NULL);
353:   ncmpi_put_vara_double_all(ncid,xyz_id,(const size_t*)&xstart,(const size_t*)&xyz_n,xarray);
354:   VecRestoreArray(xyz,&xarray);
355:   /* destroy the vectors and da */
356:   VecDestroy(natural);
357:   VecDestroy(xyz);
358:   DADestroy(dac);

360:   return(0);
361: #else /* !defined(PETSC_HAVE_PNETCDF) */
363:   SETERRQ(PETSC_ERR_SUP_SYS,"Build PETSc with NETCDF to use this viewer");
364: #endif    
365: }

367: EXTERN PetscErrorCode VecView_MPI_Draw_DA1d(Vec,PetscViewer);

372: PetscErrorCode VecView_MPI_DA(Vec xin,PetscViewer viewer)
373: {
374:   DA         da;
376:   PetscInt dim;
377:   Vec        natural;
378:   PetscTruth isdraw,ishdf4,isnetcdf;
379:   char       *prefix;

382:   PetscObjectQuery((PetscObject)xin,"DA",(PetscObject*)&da);
383:   if (!da) SETERRQ(PETSC_ERR_ARG_WRONG,"Vector not generated from a DA");
384:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
385:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_HDF4,&ishdf4);
386:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_NETCDF,&isnetcdf);
387:   if (isdraw) {
388:     DAGetInfo(da,&dim,0,0,0,0,0,0,0,0,0,0);
389:     if (dim == 1) {
390:       VecView_MPI_Draw_DA1d(xin,viewer);
391:     } else if (dim == 2) {
392:       VecView_MPI_Draw_DA2d(xin,viewer);
393:     } else {
394:       SETERRQ1(PETSC_ERR_SUP,"Cannot graphically view vector associated with this dimensional DA %D",dim);
395:     }
396:   } else if (ishdf4) {
397:     DAGetInfo(da,&dim,0,0,0,0,0,0,0,0,0,0);
398:     switch (dim) {
399:     case 2:
400:       VecView_MPI_HDF4_DA2d(xin,viewer);
401:       break;
402:     default:
403:       SETERRQ1(PETSC_ERR_SUP,"Cannot view HDF4 vector associated with this dimensional DA %D",dim);
404:     }
405:   } else if (isnetcdf) {
406:     VecView_MPI_Netcdf_DA(xin,viewer);
407:   } else {
408:     /* call viewer on natural ordering */
409:     PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);
410:     DACreateNaturalVector(da,&natural);
411:     PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);
412:     DAGlobalToNaturalBegin(da,xin,INSERT_VALUES,natural);
413:     DAGlobalToNaturalEnd(da,xin,INSERT_VALUES,natural);
414:     PetscObjectName((PetscObject)xin);
415:     PetscObjectSetName((PetscObject)natural,xin->name);
416:     VecView(natural,viewer);
417:     VecDestroy(natural);
418:   }
419:   return(0);
420: }

426: PetscErrorCode VecLoadIntoVector_Binary_DA(PetscViewer viewer,Vec xin)
427: {
428:   DA   da;
430:   Vec  natural;
431:   char *prefix;

434:   PetscObjectQuery((PetscObject)xin,"DA",(PetscObject*)&da);
435:   if (!da) SETERRQ(PETSC_ERR_ARG_WRONG,"Vector not generated from a DA");
436:   PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);
437:   DACreateNaturalVector(da,&natural);
438:   PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);
439:   VecLoadIntoVector(viewer,natural);
440:   DANaturalToGlobalBegin(da,natural,INSERT_VALUES,xin);
441:   DANaturalToGlobalEnd(da,natural,INSERT_VALUES,xin);
442:   VecDestroy(natural);
443:   PetscLogInfo(xin,"VecLoadIntoVector_Binary_DA:Loading vector from natural ordering into DA\n");
444:   return(0);
445: }