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';
 16:     buf[1] = 0;
 17:   } else {
 18:     sprintf(buf, "%0.1e", (double)val);
 19:     PetscStripZerosPlus(buf);
 20:     PetscStripe0(buf);
 21:     PetscStripInitialZero(buf);
 22:     PetscStripAllZeros(buf);
 23:     PetscStripTrailingZeros(buf);
 24:   }
 25:   *p = buf;
 26:   return 0;
 27: }

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

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

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

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

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

 55: #define EPS 1.e-6

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

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

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

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

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

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

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

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

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