xref: /petsc/src/sys/classes/draw/utils/axis.c (revision e611a964e9853b74d61a56642fe9d06a6e51780f)
1 
2 #include <../src/sys/classes/draw/utils/axisimpl.h>  /*I   "petscdraw.h"  I*/
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "PetscADefLabel"
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];
14   PetscErrorCode ierr;
15 
16   PetscFunctionBegin;
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     ierr = PetscStripZerosPlus(buf);CHKERRQ(ierr);
23     ierr = PetscStripe0(buf);CHKERRQ(ierr);
24     ierr = PetscStripInitialZero(buf);CHKERRQ(ierr);
25     ierr = PetscStripAllZeros(buf);CHKERRQ(ierr);
26     ierr = PetscStripTrailingZeros(buf);CHKERRQ(ierr);
27   }
28   *p = buf;
29   PetscFunctionReturn(0);
30 }
31 
32 #undef __FUNCT__
33 #define __FUNCT__ "PetscADefTicks"
34 /* Finds "nice" locations for the ticks */
35 PetscErrorCode PetscADefTicks(PetscReal low,PetscReal high,int num,int *ntick,PetscReal *tickloc,int maxtick)
36 {
37   PetscErrorCode ierr;
38   int            i,power;
39   PetscReal      x = 0.0,base=0.0,eps;
40 
41   PetscFunctionBegin;
42   ierr = PetscAGetBase(low,high,num,&base,&power);CHKERRQ(ierr);
43   ierr = PetscAGetNice(low,base,-1,&x);CHKERRQ(ierr);
44 
45   /* Values are of the form j * base */
46   /* Find the starting value */
47   if (x < low) x += base;
48 
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);
55 
56   if (i < 2 && num < 10) {
57     ierr = PetscADefTicks(low,high,num+1,ntick,tickloc,maxtick);CHKERRQ(ierr);
58   }
59   PetscFunctionReturn(0);
60 }
61 
62 #define EPS 1.e-6
63 
64 #undef __FUNCT__
65 #define __FUNCT__ "PetscExp10"
66 PetscErrorCode PetscExp10(PetscReal d,PetscReal *result)
67 {
68   PetscFunctionBegin;
69   *result = PetscPowReal((PetscReal)10.0,d);
70   PetscFunctionReturn(0);
71 }
72 
73 #undef __FUNCT__
74 #define __FUNCT__ "PetscMod"
75 PetscErrorCode PetscMod(PetscReal x,PetscReal y,PetscReal *result)
76 {
77   int i;
78 
79   PetscFunctionBegin;
80   if (y == 1) {
81     *result = 0.0;
82     PetscFunctionReturn(0);
83   }
84   i = ((int)x) / ((int)y);
85   x = x - i * y;
86   while (x > y) x -= y;
87   *result = x;
88   PetscFunctionReturn(0);
89 }
90 
91 #undef __FUNCT__
92 #define __FUNCT__ "PetscCopysign"
93 PetscErrorCode PetscCopysign(PetscReal a,PetscReal b,PetscReal *result)
94 {
95   PetscFunctionBegin;
96   if (b >= 0) *result = a;
97   else        *result = -a;
98   PetscFunctionReturn(0);
99 }
100 
101 #undef __FUNCT__
102 #define __FUNCT__ "PetscAGetNice"
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;
110   PetscErrorCode ierr;
111 
112   PetscFunctionBegin;
113   ierr    = PetscCopysign (0.5,(double)sign,&s);CHKERRQ(ierr);
114   etmp    = in / base + 0.5 + s;
115   ierr    = PetscCopysign (0.5,etmp,&s);CHKERRQ(ierr);
116   ierr    = PetscCopysign (EPS * etmp,(double)sign,&s2);CHKERRQ(ierr);
117   etmp    = etmp - 0.5 + s - s2;
118   ierr    = PetscMod(etmp,1.0,&m);CHKERRQ(ierr);
119   etmp    = base * (etmp -  m);
120   *result = etmp;
121   PetscFunctionReturn(0);
122 }
123 
124 #undef __FUNCT__
125 #define __FUNCT__ "PetscAGetBase"
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;
132 
133   PetscFunctionBegin;
134   /* labels of the form n * BASE */
135   /* get an approximate value for BASE */
136   base = (vmax - vmin) / (double)(num + 1);
137 
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   ierr   = PetscExp10((double)-*power,&e10);CHKERRQ(ierr);
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       ierr = PetscExp10((double)*power,&e10);CHKERRQ(ierr);
153       base = base_try[i-1] * e10;
154       if (i == 1) *power = *power + 1;
155       break;
156     }
157   }
158   *Base = base;
159   PetscFunctionReturn(0);
160 }
161