xref: /petsc/src/sys/classes/draw/utils/axis.c (revision d8e47b638cf8f604a99e9678e1df24f82d959cd7)
1 #include <petsc/private/drawimpl.h> /*I   "petscdraw.h"  I*/
2 
3 /*
4    val is the label value.  sep is the separation to the next (or previous)
5    label; this is useful in determining how many significant figures to
6    keep.
7  */
PetscADefLabel(PetscReal val,PetscReal sep,char ** p)8 PetscErrorCode PetscADefLabel(PetscReal val, PetscReal sep, char **p)
9 {
10   static char buf[40];
11 
12   PetscFunctionBegin;
13   /* Find the string */
14   if (PetscAbsReal(val) / sep < 1.e-4) {
15     buf[0] = '0';
16     buf[1] = 0;
17   } else {
18     PetscCall(PetscSNPrintf(buf, PETSC_STATIC_ARRAY_LENGTH(buf), "%0.1e", (double)val));
19     PetscCall(PetscStripZerosPlus(buf));
20     PetscCall(PetscStripe0(buf));
21     PetscCall(PetscStripInitialZero(buf));
22     PetscCall(PetscStripAllZeros(buf));
23     PetscCall(PetscStripTrailingZeros(buf));
24   }
25   *p = buf;
26   PetscFunctionReturn(PETSC_SUCCESS);
27 }
28 
29 /* Finds "nice" locations for the ticks */
PetscADefTicks(PetscReal low,PetscReal high,int num,int * ntick,PetscReal * tickloc,int maxtick)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;
34 
35   PetscFunctionBegin;
36   PetscCall(PetscAGetBase(low, high, num, &base, &power));
37   PetscCall(PetscAGetNice(low, base, -1, &x));
38 
39   /* Values are of the form j * base */
40   /* Find the starting value */
41   if (x < low) x += base;
42 
43   i   = 0;
44   eps = base / 10;
45   while (i < maxtick && x <= high + eps) {
46     tickloc[i++] = x;
47     x += base;
48   }
49   *ntick         = i;
50   tickloc[i - 1] = PetscMin(tickloc[i - 1], high);
51 
52   if (i < 2 && num < 10) PetscCall(PetscADefTicks(low, high, num + 1, ntick, tickloc, maxtick));
53   PetscFunctionReturn(PETSC_SUCCESS);
54 }
55 
56 #define EPS 1.e-6
57 
PetscExp10(PetscReal d,PetscReal * result)58 static PetscErrorCode PetscExp10(PetscReal d, PetscReal *result)
59 {
60   PetscFunctionBegin;
61   *result = PetscPowReal((PetscReal)10.0, d);
62   PetscFunctionReturn(PETSC_SUCCESS);
63 }
64 
PetscMod(PetscReal x,PetscReal y,PetscReal * result)65 static PetscErrorCode PetscMod(PetscReal x, PetscReal y, PetscReal *result)
66 {
67   int i;
68 
69   PetscFunctionBegin;
70   if (y == 1) {
71     *result = 0.0;
72     PetscFunctionReturn(PETSC_SUCCESS);
73   }
74   i = ((int)x) / ((int)y);
75   x = x - ((PetscReal)i) * y;
76   while (x > y) x -= y;
77   *result = x;
78   PetscFunctionReturn(PETSC_SUCCESS);
79 }
80 
PetscCopysign(PetscReal a,PetscReal b,PetscReal * result)81 static PetscErrorCode PetscCopysign(PetscReal a, PetscReal b, PetscReal *result)
82 {
83   PetscFunctionBegin;
84   if (b >= 0) *result = a;
85   else *result = -a;
86   PetscFunctionReturn(PETSC_SUCCESS);
87 }
88 
89 /*
90     Given a value "in" and a "base", return a nice value.
91     based on "sign", extend up (+1) or down (-1)
92  */
PetscAGetNice(PetscReal in,PetscReal base,int sign,PetscReal * result)93 PetscErrorCode PetscAGetNice(PetscReal in, PetscReal base, int sign, PetscReal *result)
94 {
95   PetscReal etmp, s, s2, m;
96 
97   PetscFunctionBegin;
98   PetscCall(PetscCopysign(0.5, (double)sign, &s));
99   etmp = in / base + 0.5 + s;
100   PetscCall(PetscCopysign(0.5, etmp, &s));
101   PetscCall(PetscCopysign(EPS * etmp, (double)sign, &s2));
102   etmp = etmp - 0.5 + s - s2;
103   PetscCall(PetscMod(etmp, 1.0, &m));
104   etmp    = base * (etmp - m);
105   *result = etmp;
106   PetscFunctionReturn(PETSC_SUCCESS);
107 }
108 
PetscAGetBase(PetscReal vmin,PetscReal vmax,int num,PetscReal * Base,int * power)109 PetscErrorCode PetscAGetBase(PetscReal vmin, PetscReal vmax, int num, PetscReal *Base, int *power)
110 {
111   PetscReal        base, ftemp, e10;
112   static PetscReal base_try[5] = {10.0, 5.0, 2.0, 1.0, 0.5};
113   int              i;
114 
115   PetscFunctionBegin;
116   /* labels of the form n * BASE */
117   /* get an approximate value for BASE */
118   base = (vmax - vmin) / (double)(num + 1);
119 
120   /* make it of form   m x 10^power,  m in [1.0, 10) */
121   if (base <= 0.0) {
122     base = PetscAbsReal(vmin);
123     if (base < 1.0) base = 1.0;
124   }
125   ftemp = PetscLog10Real((1.0 + EPS) * base);
126   if (ftemp < 0.0) ftemp -= 1.0;
127   *power = (int)ftemp;
128   PetscCall(PetscExp10((double)-*power, &e10));
129   base = base * e10;
130   if (base < 1.0) base = 1.0;
131   /* now reduce it to one of 1, 2, or 5 */
132   for (i = 1; i < 5; i++) {
133     if (base >= base_try[i]) {
134       PetscCall(PetscExp10((double)*power, &e10));
135       base = base_try[i - 1] * e10;
136       if (i == 1) *power = *power + 1;
137       break;
138     }
139   }
140   *Base = base;
141   PetscFunctionReturn(PETSC_SUCCESS);
142 }
143