Actual source code: axis.c


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

  4: /*
  5:    val is the label value.  sep is the separation to the next (or previous)
  6:    label; this is useful in determining how many significant figures to
  7:    keep.
  8:  */
  9: PetscErrorCode PetscADefLabel(PetscReal val,PetscReal sep,char **p)
 10: {
 11:   static char    buf[40];

 13:   /* Find the string */
 14:   if (PetscAbsReal(val)/sep <  1.e-4) {
 15:     buf[0] = '0'; buf[1] = 0;
 16:   } else {
 17:     sprintf(buf,"%0.1e",(double)val);
 18:     PetscStripZerosPlus(buf);
 19:     PetscStripe0(buf);
 20:     PetscStripInitialZero(buf);
 21:     PetscStripAllZeros(buf);
 22:     PetscStripTrailingZeros(buf);
 23:   }
 24:   *p = buf;
 25:   return 0;
 26: }

 28: /* Finds "nice" locations for the ticks */
 29: PetscErrorCode PetscADefTicks(PetscReal low,PetscReal high,int num,int *ntick,PetscReal *tickloc,int maxtick)
 30: {
 31:   int            i,power;
 32:   PetscReal      x = 0.0,base=0.0,eps;

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

 37:   /* Values are of the form j * base */
 38:   /* Find the starting value */
 39:   if (x < low) x += base;

 41:   i = 0; eps = base/10;
 42:   while (i < maxtick && x <= high+eps) {
 43:     tickloc[i++] = x; x += base;
 44:   }
 45:   *ntick = i;
 46:   tickloc[i-1] = PetscMin(tickloc[i-1],high);

 48:   if (i < 2 && num < 10) {
 49:     PetscADefTicks(low,high,num+1,ntick,tickloc,maxtick);
 50:   }
 51:   return 0;
 52: }

 54: #define EPS 1.e-6

 56: PetscErrorCode PetscExp10(PetscReal d,PetscReal *result)
 57: {
 58:   *result = PetscPowReal((PetscReal)10.0,d);
 59:   return 0;
 60: }

 62: PetscErrorCode PetscMod(PetscReal x,PetscReal y,PetscReal *result)
 63: {
 64:   int i;

 66:   if (y == 1) {
 67:     *result = 0.0;
 68:     return 0;
 69:   }
 70:   i = ((int)x) / ((int)y);
 71:   x = x - i * y;
 72:   while (x > y) x -= y;
 73:   *result = x;
 74:   return 0;
 75: }

 77: PetscErrorCode PetscCopysign(PetscReal a,PetscReal b,PetscReal *result)
 78: {
 79:   if (b >= 0) *result = a;
 80:   else        *result = -a;
 81:   return 0;
 82: }

 84: /*
 85:     Given a value "in" and a "base", return a nice value.
 86:     based on "sign", extend up (+1) or down (-1)
 87:  */
 88: PetscErrorCode PetscAGetNice(PetscReal in,PetscReal base,int sign,PetscReal *result)
 89: {
 90:   PetscReal      etmp,s,s2,m;

 92:   PetscCopysign (0.5,(double)sign,&s);
 93:   etmp    = in / base + 0.5 + s;
 94:   PetscCopysign (0.5,etmp,&s);
 95:   PetscCopysign (EPS * etmp,(double)sign,&s2);
 96:   etmp    = etmp - 0.5 + s - s2;
 97:   PetscMod(etmp,1.0,&m);
 98:   etmp    = base * (etmp -  m);
 99:   *result = etmp;
100:   return 0;
101: }

103: PetscErrorCode PetscAGetBase(PetscReal vmin,PetscReal vmax,int num,PetscReal *Base,int *power)
104: {
105:   PetscReal        base,ftemp,e10;
106:   static PetscReal base_try[5] = {10.0,5.0,2.0,1.0,0.5};
107:   int              i;

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

113:   /* make it of form   m x 10^power,  m in [1.0, 10) */
114:   if (base <= 0.0) {
115:     base = PetscAbsReal(vmin);
116:     if (base < 1.0) base = 1.0;
117:   }
118:   ftemp = PetscLog10Real((1.0 + EPS) * base);
119:   if (ftemp < 0.0) ftemp -= 1.0;
120:   *power = (int)ftemp;
121:   PetscExp10((double)-*power,&e10);
122:   base   = base * e10;
123:   if (base < 1.0) base = 1.0;
124:   /* now reduce it to one of 1, 2, or 5 */
125:   for (i=1; i<5; i++) {
126:     if (base >= base_try[i]) {
127:       PetscExp10((double)*power,&e10);
128:       base = base_try[i-1] * e10;
129:       if (i == 1) *power = *power + 1;
130:       break;
131:     }
132:   }
133:   *Base = base;
134:   return 0;
135: }