Actual source code: axis.c

petsc-3.3-p7 2013-05-11
  2: #include <../src/sys/draw/utils/axisimpl.h>

  6: static PetscErrorCode PetscRint(PetscReal x,PetscReal *result)
  7: {
  9:   if (x > 0) *result = floor(x + 0.5);
 10:   else       *result = floor(x - 0.5);
 11:   return(0);
 12: }

 16: /*@
 17:     PetscDrawAxisSetLimits -  Sets the limits (in user coords) of the axis
 18:     
 19:     Not Collective (ignored on all processors except processor 0 of PetscDrawAxis)

 21:     Input Parameters:
 22: +   axis - the axis
 23: .   xmin,xmax - limits in x
 24: -   ymin,ymax - limits in y

 26:     Level: advanced

 28: .seealso:  PetscDrawAxisSetHoldLimits()

 30: @*/
 31: PetscErrorCode  PetscDrawAxisSetLimits(PetscDrawAxis axis,PetscReal xmin,PetscReal xmax,PetscReal ymin,PetscReal ymax)
 32: {
 34:   if (!axis) return(0);
 35:   if (axis->hold) return(0);
 36:   axis->xlow = xmin;
 37:   axis->xhigh= xmax;
 38:   axis->ylow = ymin;
 39:   axis->yhigh= ymax;
 40:   return(0);
 41: }

 45: /*
 46:    val is the label value.  sep is the separation to the next (or previous)
 47:    label; this is useful in determining how many significant figures to   
 48:    keep.
 49:  */
 50: PetscErrorCode PetscADefLabel(PetscReal val,PetscReal sep,char **p)
 51: {
 52:   static char    buf[40];
 53:   char           fmat[10];
 55:   int            w,d;
 56:   PetscReal      rval;

 59:   /* Find the string */
 60:   if (PetscAbsReal(val)/sep <  1.e-6) {
 61:     buf[0] = '0'; buf[1] = 0;
 62:   } else if (PetscAbsReal(val) < 1.0e6 && PetscAbsReal(val) > 1.e-4) {
 63:     /* Compute the number of digits */
 64:     w = 0;
 65:     d = 0;
 66:     if (sep > 0.0) {
 67:         d = (int)ceil(- log10 (sep));
 68:         if (d < 0) d = 0;
 69:         if (PetscAbsReal(val) < 1.0e-6*sep) {
 70:             /* This is the case where we are near zero and less than a small
 71:                fraction of the sep.  In this case, we use 0 as the value */
 72:             val = 0.0;
 73:             w   = d;
 74:         }
 75:         else if (val == 0.0) w   = d;
 76:         else w = (int)(ceil(log10(PetscAbsReal(val))) + d);
 77:         if (w < 1)   w ++;
 78:         if (val < 0) w ++;
 79:     }

 81:     PetscRint(val,&rval);
 82:     if (rval == val) {
 83:         if (w > 0) sprintf(fmat,"%%%dd",w);
 84:         else {PetscStrcpy(fmat,"%d");}
 85:         sprintf(buf,fmat,(int)val);
 86:         PetscStripInitialZero(buf);
 87:         PetscStripAllZeros(buf);
 88:         PetscStripTrailingZeros(buf);
 89:     } else {
 90:         /* The code used here is inappropriate for a val of 0, which
 91:            tends to print with an excessive numer of digits.  In this
 92:            case, we should look at the next/previous values and 
 93:            use those widths */
 94:         if (w > 0) sprintf(fmat,"%%%d.%dlf",w + 1,d);
 95:         else {PetscStrcpy(fmat,"%lf");}
 96:         sprintf(buf,fmat,(double)val);
 97:         PetscStripInitialZero(buf);
 98:         PetscStripAllZeros(buf);
 99:         PetscStripTrailingZeros(buf);
100:     }
101:   } else {
102:     PetscSNPrintf(buf,40,"%g",(double)val);
103:     /* remove the extraneous 0 before the e */
104:     PetscStripZeros(buf);
105:     PetscStripZerosPlus(buf);
106:     PetscStripInitialZero(buf);
107:     PetscStripAllZeros(buf);
108:     PetscStripTrailingZeros(buf);
109:   }
110:   *p =buf;
111:   return(0);
112: }

116: /* Finds "nice" locations for the ticks */
117: PetscErrorCode PetscADefTicks(PetscReal low,PetscReal high,int num,int *ntick,PetscReal * tickloc,int  maxtick)
118: {
120:   int            i,power;
121:   PetscReal      x = 0.0,base=0.0;

124:   /* patch if low == high */
125:   if (low == high) {
126:     low  -= .01;
127:     high += .01;
128:   }

130:   /*  if (PetscAbsReal(low-high) < 1.e-8) {
131:     low  -= .01;
132:     high += .01;
133:   } */

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

138:   /* Values are of the form j * base */
139:   /* Find the starting value */
140:   if (x < low) x += base;

142:   i = 0;
143:   while (i < maxtick && x <= high) {
144:     tickloc[i++] = x;
145:     x += base;
146:   }
147:   *ntick = i;

149:   if (i < 2 && num < 10) {
150:     PetscADefTicks(low,high,num+1,ntick,tickloc,maxtick);
151:   }
152:   return(0);
153: }

155: #define EPS 1.e-6

159: PetscErrorCode PetscExp10(PetscReal d,PetscReal *result)
160: {
162:   *result = pow((PetscReal)10.0,d);
163:   return(0);
164: }

168: PetscErrorCode PetscMod(PetscReal x,PetscReal y,PetscReal *result)
169: {
170:   int     i;

173:   if (y == 1) {
174:     *result = 0.0;
175:     return(0);
176:   }
177:   i   = ((int)x) / ((int)y);
178:   x   = x - i * y;
179:   while (x > y) x -= y;
180:   *result = x;
181:   return(0);
182: }

186: PetscErrorCode PetscCopysign(PetscReal a,PetscReal b,PetscReal *result)
187: {
189:   if (b >= 0) *result = a;
190:   else        *result = -a;
191:   return(0);
192: }

196: /*
197:     Given a value "in" and a "base", return a nice value.
198:     based on "sign", extend up (+1) or down (-1)
199:  */
200: PetscErrorCode PetscAGetNice(PetscReal in,PetscReal base,int sign,PetscReal *result)
201: {
202:   PetscReal      etmp,s,s2,m;

206:   PetscCopysign (0.5,(double)sign,&s);
207:   etmp    = in / base + 0.5 + s;
208:   PetscCopysign (0.5,etmp,&s);
209:   PetscCopysign (EPS * etmp,(double)sign,&s2);
210:   etmp    = etmp - 0.5 + s - s2;
211:   PetscMod(etmp,1.0,&m);
212:   etmp    = base * (etmp -  m);
213:   *result = etmp;
214:   return(0);
215: }

219: PetscErrorCode PetscAGetBase(PetscReal vmin,PetscReal vmax,int num,PetscReal*Base,int*power)
220: {
221:   PetscReal        base,ftemp,e10;
222:   static PetscReal base_try[5] = {10.0,5.0,2.0,1.0,0.5};
223:   PetscErrorCode   ierr;
224:   int              i;

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

231:   /* make it of form   m x 10^power,  m in [1.0, 10) */
232:   if (base <= 0.0) {
233:     base    = PetscAbsReal(vmin);
234:     if (base < 1.0) base = 1.0;
235:   }
236:   ftemp   = log10((1.0 + EPS) * base);
237:   if (ftemp < 0.0)  ftemp   -= 1.0;
238:   *power  = (int)ftemp;
239:   PetscExp10((double)- *power,&e10);
240:   base    = base * e10;
241:   if (base < 1.0) base    = 1.0;
242:   /* now reduce it to one of 1, 2, or 5 */
243:   for (i=1; i<5; i++) {
244:     if (base >= base_try[i]) {
245:       PetscExp10((double)*power,&e10);
246:       base = base_try[i-1] * e10;
247:       if (i == 1) *power    = *power + 1;
248:       break;
249:     }
250:   }
251:   *Base   = base;
252:   return(0);
253: }