Actual source code: axisc.c

  1: #include <petsc/private/drawimpl.h>

  3: #define PETSC_DRAW_AXIS_MAX_SEGMENTS 20
  4: PetscClassId PETSC_DRAWAXIS_CLASSID = 0;

  6: /*@
  7:    PetscDrawAxisCreate - Generate the axis data structure.

  9:    Collective on PetscDraw

 11:    Input Parameters:
 12: .  win - PetscDraw object where axis to to be made

 14:    Output Parameter:
 15: .  axis - the axis datastructure

 17:    Notes:
 18:     the MPI communicator that owns the underlying draw object owns the PetscDrawAxis object, but calls to set PetscDrawAxis options are ignored by all processes
 19:           except the first MPI process in the communicator

 21:    Level: advanced

 23: .seealso: PetscDrawLGCreate(), PetscDrawLG, PetscDrawSPCreate(), PetscDrawSP, PetscDrawHGCreate(), PetscDrawHG, PetscDrawBarCreate(), PetscDrawBar, PetscDrawLGGetAxis(), PetscDrawSPGetAxis(),
 24:           PetscDrawHGGetAxis(), PetscDrawBarGetAxis(), PetscDrawAxis, PetscDrawAxisDestroy(), PetscDrawAxisSetColors(), PetscDrawAxisSetLabels(), PetscDrawAxisSetLimits(), PetscDrawAxisGetLimits(), PetscDrawAxisSetHoldLimits(),
 25:           PetscDrawAxisDraw()
 26: @*/
 27: PetscErrorCode  PetscDrawAxisCreate(PetscDraw draw,PetscDrawAxis *axis)
 28: {
 29:   PetscDrawAxis  ad;


 34:   PetscHeaderCreate(ad,PETSC_DRAWAXIS_CLASSID,"DrawAxis","Draw Axis","Draw",PetscObjectComm((PetscObject)draw),PetscDrawAxisDestroy,NULL);
 35:   PetscLogObjectParent((PetscObject)draw,(PetscObject)ad);

 37:   PetscObjectReference((PetscObject)draw);
 38:   ad->win = draw;

 40:   ad->xticks    = PetscADefTicks;
 41:   ad->yticks    = PetscADefTicks;
 42:   ad->xlabelstr = PetscADefLabel;
 43:   ad->ylabelstr = PetscADefLabel;
 44:   ad->ac        = PETSC_DRAW_BLACK;
 45:   ad->tc        = PETSC_DRAW_BLACK;
 46:   ad->cc        = PETSC_DRAW_BLACK;
 47:   ad->xlabel    = NULL;
 48:   ad->ylabel    = NULL;
 49:   ad->toplabel  = NULL;

 51:   *axis = ad;
 52:   return 0;
 53: }

 55: /*@
 56:     PetscDrawAxisDestroy - Frees the space used by an axis structure.

 58:     Collective on PetscDrawAxis

 60:     Input Parameters:
 61: .   axis - the axis context

 63:     Level: advanced

 65: .seealso: PetscDrawAxisCreate(), PetscDrawAxis
 66: @*/
 67: PetscErrorCode  PetscDrawAxisDestroy(PetscDrawAxis *axis)
 68: {
 69:   if (!*axis) return 0;
 71:   if (--((PetscObject)(*axis))->refct > 0) {*axis = NULL; return 0;}

 73:   PetscFree((*axis)->toplabel);
 74:   PetscFree((*axis)->xlabel);
 75:   PetscFree((*axis)->ylabel);
 76:   PetscDrawDestroy(&(*axis)->win);
 77:   PetscHeaderDestroy(axis);
 78:   return 0;
 79: }

 81: /*@
 82:     PetscDrawAxisSetColors -  Sets the colors to be used for the axis,
 83:                          tickmarks, and text.

 85:     Logically Collective on PetscDrawAxis

 87:     Input Parameters:
 88: +   axis - the axis
 89: .   ac - the color of the axis lines
 90: .   tc - the color of the tick marks
 91: -   cc - the color of the text strings

 93:     Level: advanced

 95: .seealso: PetscDrawAxisCreate(), PetscDrawAxis, PetscDrawAxisSetLabels(), PetscDrawAxisDraw(), PetscDrawAxisSetLimits()
 96: @*/
 97: PetscErrorCode  PetscDrawAxisSetColors(PetscDrawAxis axis,int ac,int tc,int cc)
 98: {
103:   axis->ac = ac; axis->tc = tc; axis->cc = cc;
104:   return 0;
105: }

107: /*@C
108:     PetscDrawAxisSetLabels -  Sets the x and y axis labels.

110:     Logically Collective on PetscDrawAxis

112:     Input Parameters:
113: +   axis - the axis
114: .   top - the label at the top of the image
115: -   xlabel,ylabel - the labes for the x and y axis

117:     Notes:
118:     Must be called before PetscDrawAxisDraw() or PetscDrawLGDraw()
119:            There should be no newlines in the arguments

121:     Level: advanced

123: .seealso: PetscDrawAxisCreate(), PetscDrawAxis, PetscDrawAxisSetColors(), PetscDrawAxisDraw(), PetscDrawAxisSetLimits()
124: @*/
125: PetscErrorCode  PetscDrawAxisSetLabels(PetscDrawAxis axis,const char top[],const char xlabel[],const char ylabel[])
126: {
128:   PetscFree(axis->xlabel);
129:   PetscFree(axis->ylabel);
130:   PetscFree(axis->toplabel);
131:   PetscStrallocpy(xlabel,&axis->xlabel);
132:   PetscStrallocpy(ylabel,&axis->ylabel);
133:   PetscStrallocpy(top,&axis->toplabel);
134:   return 0;
135: }

137: /*@
138:     PetscDrawAxisSetLimits -  Sets the limits (in user coords) of the axis

140:     Logically Collective on PetscDrawAxis

142:     Input Parameters:
143: +   axis - the axis
144: .   xmin,xmax - limits in x
145: -   ymin,ymax - limits in y

147:     Options Database:
148: .   -drawaxis_hold - hold the initial set of axis limits for future plotting

150:     Level: advanced

152: .seealso:  PetscDrawAxisSetHoldLimits(), PetscDrawAxisGetLimits(), PetscDrawAxisSetLabels(), PetscDrawAxisSetColors()

154: @*/
155: PetscErrorCode  PetscDrawAxisSetLimits(PetscDrawAxis axis,PetscReal xmin,PetscReal xmax,PetscReal ymin,PetscReal ymax)
156: {
158:   if (axis->hold) return 0;
159:   axis->xlow = xmin;
160:   axis->xhigh= xmax;
161:   axis->ylow = ymin;
162:   axis->yhigh= ymax;
163:   PetscOptionsHasName(((PetscObject)axis)->options,((PetscObject)axis)->prefix,"-drawaxis_hold",&axis->hold);
164:   return 0;
165: }

167: /*@
168:     PetscDrawAxisGetLimits -  Gets the limits (in user coords) of the axis

170:     Not Collective

172:     Input Parameters:
173: +   axis - the axis
174: .   xmin,xmax - limits in x
175: -   ymin,ymax - limits in y

177:     Level: advanced

179: .seealso:  PetscDrawAxisCreate(), PetscDrawAxis, PetscDrawAxisSetHoldLimits(), PetscDrawAxisSetLimits(), PetscDrawAxisSetLabels(), PetscDrawAxisSetColors()

181: @*/
182: PetscErrorCode  PetscDrawAxisGetLimits(PetscDrawAxis axis,PetscReal *xmin,PetscReal *xmax,PetscReal *ymin,PetscReal *ymax)
183: {
185:   if (xmin) *xmin = axis->xlow;
186:   if (xmax) *xmax = axis->xhigh;
187:   if (ymin) *ymin = axis->ylow;
188:   if (ymax) *ymax = axis->yhigh;
189:   return 0;
190: }

192: /*@
193:     PetscDrawAxisSetHoldLimits -  Causes an axis to keep the same limits until this is called
194:         again

196:     Logically Collective on PetscDrawAxis

198:     Input Parameters:
199: +   axis - the axis
200: -   hold - PETSC_TRUE - hold current limits, PETSC_FALSE allow limits to be changed

202:     Level: advanced

204:     Notes:
205:         Once this has been called with PETSC_TRUE the limits will not change if you call
206:      PetscDrawAxisSetLimits() until you call this with PETSC_FALSE

208: .seealso:  PetscDrawAxisCreate(), PetscDrawAxis, PetscDrawAxisGetLimits(), PetscDrawAxisSetLimits(), PetscDrawAxisSetLabels(), PetscDrawAxisSetColors()

210: @*/
211: PetscErrorCode  PetscDrawAxisSetHoldLimits(PetscDrawAxis axis,PetscBool hold)
212: {
215:   axis->hold = hold;
216:   return 0;
217: }

219: /*@
220:     PetscDrawAxisDraw - PetscDraws an axis.

222:     Collective on PetscDrawAxis

224:     Input Parameter:
225: .   axis - Axis structure

227:     Level: advanced

229:     Note:
230:     This draws the actual axis.  The limits etc have already been set.
231:     By picking special routines for the ticks and labels, special
232:     effects may be generated.  These routines are part of the Axis
233:     structure (axis).

235: .seealso:  PetscDrawAxisCreate(), PetscDrawAxis, PetscDrawAxisGetLimits(), PetscDrawAxisSetLimits(), PetscDrawAxisSetLabels(), PetscDrawAxisSetColors()

237: @*/
238: PetscErrorCode  PetscDrawAxisDraw(PetscDrawAxis axis)
239: {
240:   int            i,ntick,numx,numy,ac,tc,cc;
241:   PetscMPIInt    rank;
242:   size_t         len,ytlen=0;
243:   PetscReal      coors[4],tickloc[PETSC_DRAW_AXIS_MAX_SEGMENTS],sep,tw,th;
244:   PetscReal      xl,xr,yl,yr,dxl=0,dyl=0,dxr=0,dyr=0;
245:   char           *p;
246:   PetscDraw      draw;
247:   PetscBool      isnull;

251:   PetscDrawIsNull(axis->win,&isnull);
252:   if (isnull) return 0;
253:   MPI_Comm_rank(PetscObjectComm((PetscObject)axis),&rank);

255:   draw = axis->win;

257:   ac = axis->ac; tc = axis->tc; cc = axis->cc;
258:   if (axis->xlow == axis->xhigh) {axis->xlow -= .5; axis->xhigh += .5;}
259:   if (axis->ylow == axis->yhigh) {axis->ylow -= .5; axis->yhigh += .5;}

261:   PetscDrawCollectiveBegin(draw);
262:   if (rank) goto finally;

264:   /* get cannonical string size */
265:   PetscDrawSetCoordinates(draw,0,0,1,1);
266:   PetscDrawStringGetSize(draw,&tw,&th);
267:   /* lower spacing */
268:   if (axis->xlabelstr) dyl += 1.5*th;
269:   if (axis->xlabel)    dyl += 1.5*th;
270:   /* left spacing */
271:   if (axis->ylabelstr) dxl += 7.5*tw;
272:   if (axis->ylabel)    dxl += 2.0*tw;
273:   /* right and top spacing */
274:   if (axis->xlabelstr) dxr = 2.5*tw;
275:   if (axis->ylabelstr) dyr = 0.5*th;
276:   if (axis->toplabel)  dyr = 1.5*th;
277:   /* extra spacing */
278:   dxl += 0.7*tw; dxr += 0.5*tw;
279:   dyl += 0.2*th; dyr += 0.2*th;
280:   /* determine coordinates */
281:   xl = (dxl*axis->xhigh + dxr*axis->xlow - axis->xlow)  / (dxl + dxr - 1);
282:   xr = (dxl*axis->xhigh + dxr*axis->xlow - axis->xhigh) / (dxl + dxr - 1);
283:   yl = (dyl*axis->yhigh + dyr*axis->ylow - axis->ylow)  / (dyl + dyr - 1);
284:   yr = (dyl*axis->yhigh + dyr*axis->ylow - axis->yhigh) / (dyl + dyr - 1);
285:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
286:   PetscDrawStringGetSize(draw,&tw,&th);

288:   /* PetscDraw the axis lines */
289:   PetscDrawLine(draw,axis->xlow,axis->ylow,axis->xhigh,axis->ylow,ac);
290:   PetscDrawLine(draw,axis->xlow,axis->ylow,axis->xlow,axis->yhigh,ac);
291:   PetscDrawLine(draw,axis->xlow,axis->yhigh,axis->xhigh,axis->yhigh,ac);
292:   PetscDrawLine(draw,axis->xhigh,axis->ylow,axis->xhigh,axis->yhigh,ac);

294:   /* PetscDraw the top label */
295:   if (axis->toplabel) {
296:     PetscReal x = (axis->xlow + axis->xhigh)/2, y = axis->yhigh + 0.5*th;
297:     PetscDrawStringCentered(draw,x,y,cc,axis->toplabel);
298:   }

300:   /* PetscDraw the X ticks and labels */
301:   if (axis->xticks) {
302:     numx = (int)(.15*(axis->xhigh-axis->xlow)/tw); numx = PetscClipInterval(numx,2,6);
303:     (*axis->xticks)(axis->xlow,axis->xhigh,numx,&ntick,tickloc,PETSC_DRAW_AXIS_MAX_SEGMENTS);
304:     /* PetscDraw in tick marks */
305:     for (i=0; i<ntick; i++) {
306:       PetscDrawLine(draw,tickloc[i],axis->ylow,tickloc[i],axis->ylow+.5*th,tc);
307:       PetscDrawLine(draw,tickloc[i],axis->yhigh,tickloc[i],axis->yhigh-.5*th,tc);
308:     }
309:     /* label ticks */
310:     if (axis->xlabelstr) {
311:       for (i=0; i<ntick; i++) {
312:         if (i < ntick - 1) sep = tickloc[i+1] - tickloc[i];
313:         else if (i > 0)    sep = tickloc[i]   - tickloc[i-1];
314:         else               sep = 0.0;
315:         (*axis->xlabelstr)(tickloc[i],sep,&p);
316:         PetscDrawStringCentered(draw,tickloc[i],axis->ylow-1.5*th,cc,p);
317:       }
318:     }
319:   }
320:   if (axis->xlabel) {
321:     PetscReal x = (axis->xlow + axis->xhigh)/2, y = axis->ylow - 1.5*th;
322:     if (axis->xlabelstr) y -= 1.5*th;
323:     PetscDrawStringCentered(draw,x,y,cc,axis->xlabel);
324:   }

326:   /* PetscDraw the Y ticks and labels */
327:   if (axis->yticks) {
328:     numy = (int)(.50*(axis->yhigh-axis->ylow)/th); numy = PetscClipInterval(numy,2,6);
329:     (*axis->yticks)(axis->ylow,axis->yhigh,numy,&ntick,tickloc,PETSC_DRAW_AXIS_MAX_SEGMENTS);
330:     /* PetscDraw in tick marks */
331:     for (i=0; i<ntick; i++) {
332:       PetscDrawLine(draw,axis->xlow,tickloc[i],axis->xlow+.5*tw,tickloc[i],tc);
333:       PetscDrawLine(draw,axis->xhigh,tickloc[i],axis->xhigh-.5*tw,tickloc[i],tc);
334:     }
335:     /* label ticks */
336:     if (axis->ylabelstr) {
337:       for (i=0; i<ntick; i++) {
338:         if (i < ntick - 1) sep = tickloc[i+1] - tickloc[i];
339:         else if (i > 0)    sep = tickloc[i]   - tickloc[i-1];
340:         else               sep = 0.0;
341:         (*axis->ylabelstr)(tickloc[i],sep,&p);
342:         PetscStrlen(p,&len)); ytlen = PetscMax(ytlen,len;
343:         PetscDrawString(draw,axis->xlow-(len+.5)*tw,tickloc[i]-.5*th,cc,p);
344:       }
345:     }
346:   }
347:   if (axis->ylabel) {
348:     PetscReal x = axis->xlow - 2.0*tw, y = (axis->ylow + axis->yhigh)/2;
349:     if (axis->ylabelstr) x -= (ytlen+.5)*tw;
350:     PetscStrlen(axis->ylabel,&len);
351:     PetscDrawStringVertical(draw,x,y+len*th/2,cc,axis->ylabel);
352:   }

354:   PetscDrawGetCoordinates(draw,&coors[0],&coors[1],&coors[2],&coors[3]);
355: finally:
356:   PetscDrawCollectiveEnd(draw);
357:   MPI_Bcast(coors,4,MPIU_REAL,0,PetscObjectComm((PetscObject)draw));
358:   PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);
359:   return 0;
360: }

362: /*
363:     Removes all zeros but one from .0000
364: */
365: PetscErrorCode PetscStripe0(char *buf)
366: {
367:   size_t         n;
368:   PetscBool      flg;
369:   char           *str;

371:   PetscStrlen(buf,&n);
372:   PetscStrendswith(buf,"e00",&flg);
373:   if (flg) buf[n-3] = 0;
374:   PetscStrstr(buf,"e0",&str);
375:   if (str) {
376:     buf[n-2] = buf[n-1];
377:     buf[n-1] = 0;
378:   }
379:   PetscStrstr(buf,"e-0",&str);
380:   if (str) {
381:     buf[n-2] = buf[n-1];
382:     buf[n-1] = 0;
383:   }
384:   return 0;
385: }

387: /*
388:     Removes all zeros but one from .0000
389: */
390: PetscErrorCode PetscStripAllZeros(char *buf)
391: {
392:   size_t         i,n;

394:   PetscStrlen(buf,&n);
395:   if (buf[0] != '.') return 0;
396:   for (i=1; i<n; i++) {
397:     if (buf[i] != '0') return 0;
398:   }
399:   buf[0] = '0';
400:   buf[1] = 0;
401:   return 0;
402: }

404: /*
405:     Removes trailing zeros
406: */
407: PetscErrorCode PetscStripTrailingZeros(char *buf)
408: {
409:   char           *found;
410:   size_t         i,n,m = PETSC_MAX_INT;

412:   /* if there is an e in string DO NOT strip trailing zeros */
413:   PetscStrchr(buf,'e',&found);
414:   if (found) return 0;

416:   PetscStrlen(buf,&n);
417:   /* locate decimal point */
418:   for (i=0; i<n; i++) {
419:     if (buf[i] == '.') {m = i; break;}
420:   }
421:   /* if not decimal point then no zeros to remove */
422:   if (m == PETSC_MAX_INT) return 0;
423:   /* start at right end of string removing 0s */
424:   for (i=n-1; i>m; i++) {
425:     if (buf[i] != '0') return 0;
426:     buf[i] = 0;
427:   }
428:   return 0;
429: }

431: /*
432:     Removes leading 0 from 0.22 or -0.22
433: */
434: PetscErrorCode PetscStripInitialZero(char *buf)
435: {
436:   size_t         i,n;

438:   PetscStrlen(buf,&n);
439:   if (buf[0] == '0') {
440:     for (i=0; i<n; i++) buf[i] = buf[i+1];
441:   } else if (buf[0] == '-' && buf[1] == '0') {
442:     for (i=1; i<n; i++) buf[i] = buf[i+1];
443:   }
444:   return 0;
445: }

447: /*
448:      Removes the extraneous zeros in numbers like 1.10000e6
449: */
450: PetscErrorCode PetscStripZeros(char *buf)
451: {
452:   size_t         i,j,n;

454:   PetscStrlen(buf,&n);
455:   if (n<5) return 0;
456:   for (i=1; i<n-1; i++) {
457:     if (buf[i] == 'e' && buf[i-1] == '0') {
458:       for (j=i; j<n+1; j++) buf[j-1] = buf[j];
459:       PetscStripZeros(buf);
460:       return 0;
461:     }
462:   }
463:   return 0;
464: }

466: /*
467:       Removes the plus in something like 1.1e+2 or 1.1e+02
468: */
469: PetscErrorCode PetscStripZerosPlus(char *buf)
470: {
471:   size_t         i,j,n;

473:   PetscStrlen(buf,&n);
474:   if (n<5) return 0;
475:   for (i=1; i<n-2; i++) {
476:     if (buf[i] == '+') {
477:       if (buf[i+1] == '0') {
478:         for (j=i+1; j<n; j++) buf[j-1] = buf[j+1];
479:         return 0;
480:       } else {
481:         for (j=i+1; j<n+1; j++) buf[j-1] = buf[j];
482:         return 0;
483:       }
484:     } else if (buf[i] == '-') {
485:       if (buf[i+1] == '0') {
486:         for (j=i+1; j<n; j++) buf[j] = buf[j+1];
487:         return 0;
488:       }
489:     }
490:   }
491:   return 0;
492: }