Actual source code: axis.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

  2:  #include <../src/sys/classes/draw/utils/axisimpl.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];

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

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

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

 41:   /* Values are of the form j * base */
 42:   /* Find the starting value */
 43:   if (x < low) x += base;

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

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

 58: #define EPS 1.e-6

 60: PetscErrorCode PetscExp10(PetscReal d,PetscReal *result)
 61: {
 63:   *result = PetscPowReal((PetscReal)10.0,d);
 64:   return(0);
 65: }

 67: PetscErrorCode PetscMod(PetscReal x,PetscReal y,PetscReal *result)
 68: {
 69:   int i;

 72:   if (y == 1) {
 73:     *result = 0.0;
 74:     return(0);
 75:   }
 76:   i = ((int)x) / ((int)y);
 77:   x = x - i * y;
 78:   while (x > y) x -= y;
 79:   *result = x;
 80:   return(0);
 81: }

 83: PetscErrorCode PetscCopysign(PetscReal a,PetscReal b,PetscReal *result)
 84: {
 86:   if (b >= 0) *result = a;
 87:   else        *result = -a;
 88:   return(0);
 89: }

 91: /*
 92:     Given a value "in" and a "base", return a nice value.
 93:     based on "sign", extend up (+1) or down (-1)
 94:  */
 95: PetscErrorCode PetscAGetNice(PetscReal in,PetscReal base,int sign,PetscReal *result)
 96: {
 97:   PetscReal      etmp,s,s2,m;

101:   PetscCopysign (0.5,(double)sign,&s);
102:   etmp    = in / base + 0.5 + s;
103:   PetscCopysign (0.5,etmp,&s);
104:   PetscCopysign (EPS * etmp,(double)sign,&s2);
105:   etmp    = etmp - 0.5 + s - s2;
106:   PetscMod(etmp,1.0,&m);
107:   etmp    = base * (etmp -  m);
108:   *result = etmp;
109:   return(0);
110: }

112: PetscErrorCode PetscAGetBase(PetscReal vmin,PetscReal vmax,int num,PetscReal *Base,int *power)
113: {
114:   PetscReal        base,ftemp,e10;
115:   static PetscReal base_try[5] = {10.0,5.0,2.0,1.0,0.5};
116:   PetscErrorCode   ierr;
117:   int              i;

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

124:   /* make it of form   m x 10^power,  m in [1.0, 10) */
125:   if (base <= 0.0) {
126:     base = PetscAbsReal(vmin);
127:     if (base < 1.0) base = 1.0;
128:   }
129:   ftemp = PetscLog10Real((1.0 + EPS) * base);
130:   if (ftemp < 0.0) ftemp -= 1.0;
131:   *power = (int)ftemp;
132:   PetscExp10((double)-*power,&e10);
133:   base   = base * e10;
134:   if (base < 1.0) base = 1.0;
135:   /* now reduce it to one of 1, 2, or 5 */
136:   for (i=1; i<5; i++) {
137:     if (base >= base_try[i]) {
138:       PetscExp10((double)*power,&e10);
139:       base = base_try[i-1] * e10;
140:       if (i == 1) *power = *power + 1;
141:       break;
142:     }
143:   }
144:   *Base = base;
145:   return(0);
146: }