Actual source code: axis.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  2: #include <../src/sys/classes/draw/utils/axisimpl.h>  /*I   "petscdraw.h"  I*/

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

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

 34: /* Finds "nice" locations for the ticks */
 35: PetscErrorCode PetscADefTicks(PetscReal low,PetscReal high,int num,int *ntick,PetscReal *tickloc,int maxtick)
 36: {
 38:   int            i,power;
 39:   PetscReal      x = 0.0,base=0.0,eps;

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

 45:   /* Values are of the form j * base */
 46:   /* Find the starting value */
 47:   if (x < low) x += base;

 49:   i = 0; eps = base/10;
 50:   while (i < maxtick && x <= high+eps) {
 51:     tickloc[i++] = x; x += base;
 52:   }
 53:   *ntick = i;
 54:   tickloc[i-1] = PetscMin(tickloc[i-1],high);

 56:   if (i < 2 && num < 10) {
 57:     PetscADefTicks(low,high,num+1,ntick,tickloc,maxtick);
 58:   }
 59:   return(0);
 60: }

 62: #define EPS 1.e-6

 66: PetscErrorCode PetscExp10(PetscReal d,PetscReal *result)
 67: {
 69:   *result = PetscPowReal((PetscReal)10.0,d);
 70:   return(0);
 71: }

 75: PetscErrorCode PetscMod(PetscReal x,PetscReal y,PetscReal *result)
 76: {
 77:   int i;

 80:   if (y == 1) {
 81:     *result = 0.0;
 82:     return(0);
 83:   }
 84:   i = ((int)x) / ((int)y);
 85:   x = x - i * y;
 86:   while (x > y) x -= y;
 87:   *result = x;
 88:   return(0);
 89: }

 93: PetscErrorCode PetscCopysign(PetscReal a,PetscReal b,PetscReal *result)
 94: {
 96:   if (b >= 0) *result = a;
 97:   else        *result = -a;
 98:   return(0);
 99: }

103: /*
104:     Given a value "in" and a "base", return a nice value.
105:     based on "sign", extend up (+1) or down (-1)
106:  */
107: PetscErrorCode PetscAGetNice(PetscReal in,PetscReal base,int sign,PetscReal *result)
108: {
109:   PetscReal      etmp,s,s2,m;

113:   PetscCopysign (0.5,(double)sign,&s);
114:   etmp    = in / base + 0.5 + s;
115:   PetscCopysign (0.5,etmp,&s);
116:   PetscCopysign (EPS * etmp,(double)sign,&s2);
117:   etmp    = etmp - 0.5 + s - s2;
118:   PetscMod(etmp,1.0,&m);
119:   etmp    = base * (etmp -  m);
120:   *result = etmp;
121:   return(0);
122: }

126: PetscErrorCode PetscAGetBase(PetscReal vmin,PetscReal vmax,int num,PetscReal *Base,int *power)
127: {
128:   PetscReal        base,ftemp,e10;
129:   static PetscReal base_try[5] = {10.0,5.0,2.0,1.0,0.5};
130:   PetscErrorCode   ierr;
131:   int              i;

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

138:   /* make it of form   m x 10^power,  m in [1.0, 10) */
139:   if (base <= 0.0) {
140:     base = PetscAbsReal(vmin);
141:     if (base < 1.0) base = 1.0;
142:   }
143:   ftemp = PetscLog10Real((1.0 + EPS) * base);
144:   if (ftemp < 0.0) ftemp -= 1.0;
145:   *power = (int)ftemp;
146:   PetscExp10((double)-*power,&e10);
147:   base   = base * e10;
148:   if (base < 1.0) base = 1.0;
149:   /* now reduce it to one of 1, 2, or 5 */
150:   for (i=1; i<5; i++) {
151:     if (base >= base_try[i]) {
152:       PetscExp10((double)*power,&e10);
153:       base = base_try[i-1] * e10;
154:       if (i == 1) *power = *power + 1;
155:       break;
156:     }
157:   }
158:   *Base = base;
159:   return(0);
160: }