Actual source code: axis.c

  1: /*$Id: axis.c,v 1.76 2001/08/07 21:28:45 bsmith Exp $*/
  2: /*
  3:    This file contains a simple routine for generating a 2-d axis.
  4: */

 6:  #include petsc.h

  8: int DRAWAXIS_COOKIE;

 10: struct _p_DrawAxis {
 11:     PETSCHEADER(int)
 12:     PetscReal  xlow,ylow,xhigh,yhigh;     /* User - coord limits */
 13:     int        (*ylabelstr)(PetscReal,PetscReal,char **),/* routines to generate labels */
 14:                (*xlabelstr)(PetscReal,PetscReal,char **);
 15:     int        (*xticks)(PetscReal,PetscReal,int,int*,PetscReal*,int),
 16:                (*yticks)(PetscReal,PetscReal,int,int*,PetscReal*,int);
 17:                                           /* location and size of ticks */
 18:     PetscDraw  win;
 19:     int        ac,tc,cc;                     /* axis,tick, character color */
 20:     char       *xlabel,*ylabel,*toplabel;
 21:     PetscTruth hold;
 22: };

 24: #define MAXSEGS 20

 26: EXTERN int    PetscADefTicks(PetscReal,PetscReal,int,int*,PetscReal*,int);
 27: EXTERN int    PetscADefLabel(PetscReal,PetscReal,char**);
 28: static int    PetscAGetNice(PetscReal,PetscReal,int,PetscReal*);
 29: static int    PetscAGetBase(PetscReal,PetscReal,int,PetscReal*,int*);

 31: #undef __FUNCT__  
 33: static int PetscRint(PetscReal x,PetscReal *result)
 34: {
 36:   if (x > 0) *result = floor(x + 0.5);
 37:   else       *result = floor(x - 0.5);
 38:   return(0);
 39: }

 41: #undef __FUNCT__  
 43: /*@C
 44:    PetscDrawAxisCreate - Generate the axis data structure.

 46:    Collective over PetscDraw

 48:    Input Parameters:
 49: .  win - PetscDraw object where axis to to be made

 51:    Ouput Parameters:
 52: .  axis - the axis datastructure

 54:    Level: advanced

 56: @*/
 57: int PetscDrawAxisCreate(PetscDraw draw,PetscDrawAxis *axis)
 58: {
 59:   PetscDrawAxis ad;
 60:   PetscObject   obj = (PetscObject)draw;
 61:   int           ierr;
 62:   PetscTruth    isnull;

 67:   PetscTypeCompare(obj,PETSC_DRAW_NULL,&isnull);
 68:   if (isnull) {
 69:     PetscDrawOpenNull(obj->comm,(PetscDraw*)axis);
 70:     (*axis)->win = draw;
 71:     return(0);
 72:   }
 73:   PetscHeaderCreate(ad,_p_DrawAxis,int,DRAWAXIS_COOKIE,0,"PetscDrawAxis",obj->comm,PetscDrawAxisDestroy,0);
 74:   PetscLogObjectCreate(ad);
 75:   PetscLogObjectParent(draw,ad);
 76:   ad->xticks    = PetscADefTicks;
 77:   ad->yticks    = PetscADefTicks;
 78:   ad->xlabelstr = PetscADefLabel;
 79:   ad->ylabelstr = PetscADefLabel;
 80:   ad->win       = draw;
 81:   ad->ac        = PETSC_DRAW_BLACK;
 82:   ad->tc        = PETSC_DRAW_BLACK;
 83:   ad->cc        = PETSC_DRAW_BLACK;
 84:   ad->xlabel    = 0;
 85:   ad->ylabel    = 0;
 86:   ad->toplabel  = 0;

 88:   *axis = ad;
 89:   return(0);
 90: }

 92: #undef __FUNCT__  
 94: /*@C
 95:     PetscDrawAxisDestroy - Frees the space used by an axis structure.

 97:     Collective over PetscDrawAxis

 99:     Input Parameters:
100: .   axis - the axis context
101:  
102:     Level: advanced

104: @*/
105: int PetscDrawAxisDestroy(PetscDrawAxis axis)
106: {

110:   if (!axis) return(0);
111:   if (--axis->refct > 0) return(0);

113:   PetscStrfree(axis->toplabel);
114:   PetscStrfree(axis->xlabel);
115:   PetscStrfree(axis->ylabel);
116:   PetscLogObjectDestroy(axis);
117:   PetscHeaderDestroy(axis);
118:   return(0);
119: }

121: #undef __FUNCT__  
123: /*@
124:     PetscDrawAxisSetColors -  Sets the colors to be used for the axis,       
125:                          tickmarks, and text.

127:     Not Collective (ignored on all processors except processor 0 of PetscDrawAxis)

129:     Input Parameters:
130: +   axis - the axis
131: .   ac - the color of the axis lines
132: .   tc - the color of the tick marks
133: -   cc - the color of the text strings

135:     Level: advanced

137: @*/
138: int PetscDrawAxisSetColors(PetscDrawAxis axis,int ac,int tc,int cc)
139: {
141:   if (!axis) return(0);
142:   axis->ac = ac; axis->tc = tc; axis->cc = cc;
143:   return(0);
144: }

146: #undef __FUNCT__  
148: /*@C
149:     PetscDrawAxisSetLabels -  Sets the x and y axis labels.

151:     Not Collective (ignored on all processors except processor 0 of PetscDrawAxis)

153:     Input Parameters:
154: +   axis - the axis
155: .   top - the label at the top of the image
156: -   xlabel,ylabel - the labes for the x and y axis

158:     Level: advanced

160: @*/
161: int PetscDrawAxisSetLabels(PetscDrawAxis axis,char* top,char *xlabel,char *ylabel)
162: {

166:   if (!axis) return(0);
167:   PetscStrallocpy(xlabel,&axis->xlabel);
168:   PetscStrallocpy(ylabel,&axis->ylabel);
169:   PetscStrallocpy(top,&axis->toplabel);
170:   return(0);
171: }

173: #undef __FUNCT__  
175: /*@
176:     PetscDrawAxisSetHoldLimits -  Causes an axis to keep the same limits until this is called
177:         again
178:     
179:     Not Collective (ignored on all processors except processor 0 of PetscDrawAxis)

181:     Input Parameters:
182: +   axis - the axis
183: -   hold - PETSC_TRUE - hold current limits, PETSC_FALSE allow limits to be changed

185:     Level: advanced

187:     Notes:
188:         Once this has been called with PETSC_TRUE the limits will not change if you call
189:      PetscDrawAxisSetLimits() until you call this with PETSC_FALSE
190:  
191: .seealso:  PetscDrawAxisSetLimits()

193: @*/
194: int PetscDrawAxisSetHoldLimits(PetscDrawAxis axis,PetscTruth hold)
195: {
197:   if (!axis) return(0);
198:   axis->hold = hold;
199:   return(0);
200: }

202: #undef __FUNCT__  
204: /*@
205:     PetscDrawAxisSetLimits -  Sets the limits (in user coords) of the axis
206:     
207:     Not Collective (ignored on all processors except processor 0 of PetscDrawAxis)

209:     Input Parameters:
210: +   axis - the axis
211: .   xmin,xmax - limits in x
212: -   ymin,ymax - limits in y

214:     Level: advanced

216: .seealso:  PetscDrawAxisSetHoldLimits()

218: @*/
219: int PetscDrawAxisSetLimits(PetscDrawAxis axis,PetscReal xmin,PetscReal xmax,PetscReal ymin,PetscReal ymax)
220: {
222:   if (!axis) return(0);
223:   if (axis->hold) return(0);
224:   axis->xlow = xmin;
225:   axis->xhigh= xmax;
226:   axis->ylow = ymin;
227:   axis->yhigh= ymax;
228:   return(0);
229: }

231: #undef __FUNCT__  
233: /*@
234:     PetscDrawAxisDraw - PetscDraws an axis.

236:     Not Collective (ignored on all processors except processor 0 of PetscDrawAxis)

238:     Input Parameter:
239: .   axis - Axis structure

241:     Level: advanced

243:     Note:
244:     This draws the actual axis.  The limits etc have already been set.
245:     By picking special routines for the ticks and labels, special
246:     effects may be generated.  These routines are part of the Axis
247:     structure (axis).
248: @*/
249: int PetscDrawAxisDraw(PetscDrawAxis axis)
250: {
251:   int       i,ierr,ntick,numx,numy,ac = axis->ac,tc = axis->tc,cc = axis->cc,rank,len;
252:   PetscReal tickloc[MAXSEGS],sep,h,w,tw,th,xl,xr,yl,yr;
253:   char      *p;
254:   PetscDraw draw = axis->win;
255: 
257:   if (!axis) return(0);
258:   MPI_Comm_rank(axis->comm,&rank);
259:   if (rank) return(0);

261:   if (axis->xlow == axis->xhigh) {axis->xlow -= .5; axis->xhigh += .5;}
262:   if (axis->ylow == axis->yhigh) {axis->ylow -= .5; axis->yhigh += .5;}
263:   xl = axis->xlow; xr = axis->xhigh; yl = axis->ylow; yr = axis->yhigh;
264:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
265:   PetscDrawStringGetSize(draw,&tw,&th);
266:   numx = (int)(.15*(xr-xl)/tw); if (numx > 6) numx = 6; if (numx< 2) numx = 2;
267:   numy = (int)(.5*(yr-yl)/th); if (numy > 6) numy = 6; if (numy< 2) numy = 2;
268:   xl -= 8*tw; xr += 2*tw; yl -= 2.5*th; yr += 2*th;
269:   if (axis->xlabel) yl -= 2*th;
270:   if (axis->ylabel) xl -= 2*tw;
271:   PetscDrawSetCoordinates(draw,xl,yl,xr,yr);
272:   PetscDrawStringGetSize(draw,&tw,&th);

274:   PetscDrawLine(draw,axis->xlow,axis->ylow,axis->xhigh,axis->ylow,ac);
275:   PetscDrawLine(draw,axis->xlow,axis->ylow,axis->xlow,axis->yhigh,ac);

277:   if (axis->toplabel) {
278:      PetscStrlen(axis->toplabel,&len);
279:     w    = xl + .5*(xr - xl) - .5*len*tw;
280:     h    = axis->yhigh;
281:     PetscDrawString(draw,w,h,cc,axis->toplabel);
282:   }

284:   /* PetscDraw the ticks and labels */
285:   if (axis->xticks) {
286:     (*axis->xticks)(axis->xlow,axis->xhigh,numx,&ntick,tickloc,MAXSEGS);
287:     /* PetscDraw in tick marks */
288:     for (i=0; i<ntick; i++) {
289:       PetscDrawLine(draw,tickloc[i],axis->ylow-.5*th,tickloc[i],axis->ylow+.5*th,tc);
290:     }
291:     /* label ticks */
292:     for (i=0; i<ntick; i++) {
293:         if (axis->xlabelstr) {
294:             if (i < ntick - 1) sep = tickloc[i+1] - tickloc[i];
295:             else if (i > 0)    sep = tickloc[i]   - tickloc[i-1];
296:             else               sep = 0.0;
297:             (*axis->xlabelstr)(tickloc[i],sep,&p);
298:             PetscStrlen(p,&len);
299:             w    = .5*len*tw;
300:             PetscDrawString(draw,tickloc[i]-w,axis->ylow-1.2*th,cc,p);
301:         }
302:     }
303:   }
304:   if (axis->xlabel) {
305:     PetscStrlen(axis->xlabel,&len);
306:     w    = xl + .5*(xr - xl) - .5*len*tw;
307:     h    = axis->ylow - 2.5*th;
308:     PetscDrawString(draw,w,h,cc,axis->xlabel);
309:   }
310:   if (axis->yticks) {
311:     (*axis->yticks)(axis->ylow,axis->yhigh,numy,&ntick,tickloc,MAXSEGS);
312:     /* PetscDraw in tick marks */
313:     for (i=0; i<ntick; i++) {
314:       PetscDrawLine(draw,axis->xlow -.5*tw,tickloc[i],axis->xlow+.5*tw,tickloc[i],tc);
315:     }
316:     /* label ticks */
317:     for (i=0; i<ntick; i++) {
318:         if (axis->ylabelstr) {
319:             if (i < ntick - 1) sep = tickloc[i+1] - tickloc[i];
320:             else if (i > 0)    sep = tickloc[i]   - tickloc[i-1];
321:             else               sep = 0.0;
322:             (*axis->xlabelstr)(tickloc[i],sep,&p);
323:             PetscStrlen(p,&len);
324:             w    = axis->xlow - len * tw - 1.2*tw;
325:             PetscDrawString(draw,w,tickloc[i]-.5*th,cc,p);
326:         }
327:     }
328:   }
329:   if (axis->ylabel) {
330:     PetscStrlen(axis->ylabel,&len);
331:     h    = yl + .5*(yr - yl) + .5*len*th;
332:     w    = xl + .5*tw;
333:     PetscDrawStringVertical(draw,w,h,cc,axis->ylabel);
334:   }
335:   return(0);
336: }

338: #undef __FUNCT__  
340: /*
341:     Removes all zeros but one from .0000 
342: */
343: static int PetscStripAllZeros(char *buf)
344: {
345:   int i,n,ierr;

348:   PetscStrlen(buf,&n);
349:   if (buf[0] != '.') return(0);
350:   for (i=1; i<n; i++) {
351:     if (buf[i] != '0') return(0);
352:   }
353:   buf[0] = '0';
354:   buf[1] = 0;
355:   return(0);
356: }

358: #undef __FUNCT__  
360: /*
361:     Removes trailing zeros
362: */
363: static int PetscStripTrailingZeros(char *buf)
364: {
365:   int  ierr,i,n,m = -1;
366:   char *found;

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

373:   PetscStrlen(buf,&n);
374:   /* locate decimal point */
375:   for (i=0; i<n; i++) {
376:     if (buf[i] == '.') {m = i; break;}
377:   }
378:   /* if not decimal point then no zeros to remove */
379:   if (m == -1) return(0);
380:   /* start at right end of string removing 0s */
381:   for (i=n-1; i>m; i++) {
382:     if (buf[i] != '0') return(0);
383:     buf[i] = 0;
384:   }
385:   return(0);
386: }

388: #undef __FUNCT__  
390: /*
391:     Removes leading 0 from 0.22 or -0.22
392: */
393: static int PetscStripInitialZero(char *buf)
394: {
395:   int i,n,ierr;

398:   PetscStrlen(buf,&n);
399:   if (buf[0] == '0') {
400:     for (i=0; i<n; i++) {
401:       buf[i] = buf[i+1];
402:     }
403:   } else if (buf[0] == '-' && buf[1] == '0') {
404:     for (i=1; i<n; i++) {
405:       buf[i] = buf[i+1];
406:     }
407:   }
408:   return(0);
409: }

411: #undef __FUNCT__  
413: /*
414:      Removes the extraneous zeros in numbers like 1.10000e6
415: */
416: static int PetscStripZeros(char *buf)
417: {
418:   int i,j,n,ierr;

421:   PetscStrlen(buf,&n);
422:   if (n<5) return(0);
423:   for (i=1; i<n-1; i++) {
424:     if (buf[i] == 'e' && buf[i-1] == '0') {
425:       for (j=i; j<n+1; j++) buf[j-1] = buf[j];
426:       PetscStripZeros(buf);
427:       return(0);
428:     }
429:   }
430:   return(0);
431: }

433: #undef __FUNCT__  
435: /*
436:       Removes the plus in something like 1.1e+2
437: */
438: static int PetscStripZerosPlus(char *buf)
439: {
440:   int i,j,n,ierr;

443:   PetscStrlen(buf,&n);
444:   if (n<5) return(0);
445:   for (i=1; i<n-2; i++) {
446:     if (buf[i] == '+') {
447:       if (buf[i+1] == '0') {
448:         for (j=i+1; j<n+1; j++) buf[j-1] = buf[j+1];
449:         return(0);
450:       } else {
451:         for (j=i+1; j<n+1; j++) buf[j] = buf[j+1];
452:         return(0);
453:       }
454:     } else if (buf[i] == '-') {
455:       if (buf[i+1] == '0') {
456:         for (j=i+1; j<n+1; j++) buf[j] = buf[j+1];
457:         return(0);
458:       }
459:     }
460:   }
461:   return(0);
462: }

464: #undef __FUNCT__  
466: /*
467:    val is the label value.  sep is the separation to the next (or previous)
468:    label; this is useful in determining how many significant figures to   
469:    keep.
470:  */
471: int PetscADefLabel(PetscReal val,PetscReal sep,char **p)
472: {
473:   static char buf[40];
474:   char        fmat[10];
475:   int         ierr,w,d;
476:   PetscReal   rval;

479:   /* Find the string */
480:   if (PetscAbsReal(val)/sep <  1.e-6) {
481:     buf[0] = '0'; buf[1] = 0;
482:   } else if (PetscAbsReal(val) < 1.0e6 && PetscAbsReal(val) > 1.e-4) {
483:     /* Compute the number of digits */
484:     w = 0;
485:     d = 0;
486:     if (sep > 0.0) {
487:         d = (int)ceil(- log10 (sep));
488:         if (d < 0) d = 0;
489:         if (PetscAbsReal(val) < 1.0e-6*sep) {
490:             /* This is the case where we are near zero and less than a small
491:                fraction of the sep.  In this case, we use 0 as the value */
492:             val = 0.0;
493:             w   = d;
494:         }
495:         else if (val == 0.0) w   = d;
496:         else w = (int)(ceil(log10(PetscAbsReal(val))) + d);
497:         if (w < 1)   w ++;
498:         if (val < 0) w ++;
499:     }

501:     PetscRint(val,&rval);
502:     if (rval == val) {
503:         if (w > 0) sprintf(fmat,"%%%dd",w);
504:         else {PetscStrcpy(fmat,"%d");}
505:         sprintf(buf,fmat,(int)val);
506:         PetscStripInitialZero(buf);
507:         PetscStripAllZeros(buf);
508:         PetscStripTrailingZeros(buf);
509:     } else {
510:         /* The code used here is inappropriate for a val of 0, which
511:            tends to print with an excessive numer of digits.  In this
512:            case, we should look at the next/previous values and 
513:            use those widths */
514:         if (w > 0) sprintf(fmat,"%%%d.%dlf",w + 1,d);
515:         else {PetscStrcpy(fmat,"%lf");}
516:         sprintf(buf,fmat,val);
517:         PetscStripInitialZero(buf);
518:         PetscStripAllZeros(buf);
519:         PetscStripTrailingZeros(buf);
520:     }
521:   } else {
522:     sprintf(buf,"%e",val);
523:     /* remove the extraneous 0 before the e */
524:     PetscStripZeros(buf);
525:     PetscStripZerosPlus(buf);
526:     PetscStripInitialZero(buf);
527:     PetscStripAllZeros(buf);
528:     PetscStripTrailingZeros(buf);
529:   }
530:   *p =buf;
531:   return(0);
532: }

534: #undef __FUNCT__  
536: /* Finds "nice" locations for the ticks */
537: int PetscADefTicks(PetscReal low,PetscReal high,int num,int *ntick,PetscReal * tickloc,int  maxtick)
538: {
539:   int       i,power,ierr;
540:   PetscReal x,base;

543:   /* patch if low == high */
544:   if (low == high) {
545:     low  -= .01;
546:     high += .01;
547:   }

549:   /*  if (PetscAbsReal(low-high) < 1.e-8) {
550:     low  -= .01;
551:     high += .01;
552:   } */

554:   PetscAGetBase(low,high,num,&base,&power);
555:   PetscAGetNice(low,base,-1,&x);

557:   /* Values are of the form j * base */
558:   /* Find the starting value */
559:   if (x < low) x += base;

561:   i = 0;
562:   while (i < maxtick && x <= high) {
563:     tickloc[i++] = x;
564:     x += base;
565:   }
566:   *ntick = i;

568:   if (i < 2 && num < 10) {
569:     PetscADefTicks(low,high,num+1,ntick,tickloc,maxtick);
570:   }
571:   return(0);
572: }

574: #define EPS 1.e-6

576: #undef __FUNCT__  
578: static int PetscExp10(PetscReal d,PetscReal *result)
579: {
581:   *result = pow(10.0,d);
582:   return(0);
583: }

585: #undef __FUNCT__  
587: static int PetscMod(PetscReal x,PetscReal y,PetscReal *result)
588: {
589:   int     i;

592:   i   = ((int)x) / ((int)y);
593:   x   = x - i * y;
594:   while (x > y) x -= y;
595:   *result = x;
596:   return(0);
597: }

599: #undef __FUNCT__  
601: static int PetscCopysign(PetscReal a,PetscReal b,PetscReal *result)
602: {
604:   if (b >= 0) *result = a;
605:   else        *result = -a;
606:   return(0);
607: }

609: #undef __FUNCT__  
611: /*
612:     Given a value "in" and a "base", return a nice value.
613:     based on "sign", extend up (+1) or down (-1)
614:  */
615: static int PetscAGetNice(PetscReal in,PetscReal base,int sign,PetscReal *result)
616: {
617:   PetscReal  etmp,s,s2,m;
618:   int     ierr;

621:   ierr    = PetscCopysign (0.5,(double)sign,&s);
622:   etmp    = in / base + 0.5 + s;
623:   ierr    = PetscCopysign (0.5,etmp,&s);
624:   ierr    = PetscCopysign (EPS * etmp,(double)sign,&s2);
625:   etmp    = etmp - 0.5 + s - s2;
626:   ierr    = PetscMod(etmp,1.0,&m);
627:   etmp    = base * (etmp -  m);
628:   *result = etmp;
629:   return(0);
630: }

632: #undef __FUNCT__  
634: static int PetscAGetBase(PetscReal vmin,PetscReal vmax,int num,PetscReal*Base,int*power)
635: {
636:   PetscReal        base,ftemp,e10;
637:   static PetscReal base_try[5] = {10.0,5.0,2.0,1.0,0.5};
638:   int              i,ierr;

641:   /* labels of the form n * BASE */
642:   /* get an approximate value for BASE */
643:   base    = (vmax - vmin) / (double)(num + 1);

645:   /* make it of form   m x 10^power,  m in [1.0, 10) */
646:   if (base <= 0.0) {
647:     base    = PetscAbsReal(vmin);
648:     if (base < 1.0) base = 1.0;
649:   }
650:   ftemp   = log10((1.0 + EPS) * base);
651:   if (ftemp < 0.0)  ftemp   -= 1.0;
652:   *power  = (int)ftemp;
653:   PetscExp10((double)- *power,&e10);
654:   base    = base * e10;
655:   if (base < 1.0) base    = 1.0;
656:   /* now reduce it to one of 1, 2, or 5 */
657:   for (i=1; i<5; i++) {
658:     if (base >= base_try[i]) {
659:       PetscExp10((double)*power,&e10);
660:       base = base_try[i-1] * e10;
661:       if (i == 1) *power    = *power + 1;
662:       break;
663:     }
664:   }
665:   *Base   = base;
666:   return(0);
667: }