xref: /petsc/src/sys/classes/draw/utils/cmap.c (revision 00d931fe9835bef04c3bcd2a9a1bf118d64cc4c2)
1 #include <petscsys.h>              /*I "petscsys.h" I*/
2 #include <petscdraw.h>
3 
4 /*
5     Set up a color map, using uniform separation in hue space.
6     Map entries are Red, Green, Blue.
7     Values are "gamma" corrected.
8  */
9 
10 /*
11    Gamma is a monitor dependent value.  The value here is an
12    approximate that gives somewhat better results than Gamma = 1.
13  */
14 static PetscReal Gamma = 2.0;
15 
16 #undef __FUNCT__
17 #define __FUNCT__ "PetscDrawUtilitySetGamma"
18 PetscErrorCode  PetscDrawUtilitySetGamma(PetscReal g)
19 {
20   PetscFunctionBegin;
21   Gamma = g;
22   PetscFunctionReturn(0);
23 }
24 
25 /*
26  * This algorithm is from Foley and van Dam, page 616
27  * given
28  *   (0:359, 0:100, 0:100).
29  *      h       l      s
30  * set
31  *   (0:255, 0:255, 0:255)
32  *      r       g      b
33  */
34 PETSC_STATIC_INLINE int PetscHlsHelper(int h,int m1,int m2)
35 {
36   while (h > 360) h = h - 360;
37   while (h < 0)   h = h + 360;
38   if (h < 60)  return m1 + (m2-m1)*h/60;
39   if (h < 180) return m2;
40   if (h < 240) return m1 + (m2-m1)*(240-h)/60;
41   return m1;
42 }
43 
44 PETSC_STATIC_INLINE void PetscHlsToRgb(int h,int l,int s,unsigned char *r,unsigned char *g,unsigned char *b)
45 {
46   int m1,m2;         /* in 0 to 100 */
47 
48   if (l <= 50) m2 = l * (100 + s) / 100 ; /* not sure of "/100" */
49   else         m2 = l + s - l*s/100;
50 
51   m1 = 2*l - m2;
52   if (!s) {
53     /* ignore h */
54     *r = 255 * l / 100;
55     *g = 255 * l / 100;
56     *b = 255 * l / 100;
57   } else {
58     *r = (255 * PetscHlsHelper(h+120,m1,m2)) / 100;
59     *g = (255 * PetscHlsHelper(h,m1,m2))     / 100;
60     *b = (255 * PetscHlsHelper(h-120,m1,m2)) / 100;
61   }
62 }
63 
64 #undef __FUNCT__
65 #define __FUNCT__ "PetscDrawCmap_Hue"
66 static PetscErrorCode PetscDrawCmap_Hue(int mapsize, unsigned char R[],unsigned char G[],unsigned char B[])
67 {
68   int            i,hue,lightness,saturation;
69   PetscReal      igamma = 1.0 / Gamma;
70 
71   PetscFunctionBegin;
72   hue        = 0;    /* in 0:359 */
73   lightness  = 50;   /* in 0:100 */
74   saturation = 100;  /* in 0:100 */
75   for (i=0; i<mapsize; i++) {
76     PetscHlsToRgb(hue,lightness,saturation,&R[i],&G[i],&B[i]);;
77     R[i] = (unsigned char)(PetscFloorReal(255.999*PetscPowReal(((PetscReal)R[i])/255,igamma)));
78     G[i] = (unsigned char)(PetscFloorReal(255.999*PetscPowReal(((PetscReal)G[i])/255,igamma)));
79     B[i] = (unsigned char)(PetscFloorReal(255.999*PetscPowReal(((PetscReal)B[i])/255,igamma)));
80     hue += (359/(mapsize-2));
81   }
82   PetscFunctionReturn(0);
83 }
84 
85 #undef __FUNCT__
86 #define __FUNCT__ "PetscDrawCmap_Gray"
87 static PetscErrorCode PetscDrawCmap_Gray(int mapsize,unsigned char R[],unsigned char G[],unsigned char B[])
88 {
89   int i;
90   PetscFunctionBegin;
91   for (i=0; i<mapsize; i++) R[i] = G[i] = B[i] = (unsigned char)((255.0*i)/(mapsize-1));
92   PetscFunctionReturn(0);
93 }
94 
95 #undef __FUNCT__
96 #define __FUNCT__ "PetscDrawCmap_Jet"
97 static PetscErrorCode PetscDrawCmap_Jet(int mapsize,unsigned char R[],unsigned char G[],unsigned char B[])
98 {
99   int          i;
100   const double knots[] =  {0, 1/8., 3/8., 5/8., 7/8., 1};
101 
102   PetscFunctionBegin;
103   for (i=0; i<mapsize; i++) {
104     double u = (double)i/(mapsize-1);
105     double m, r=0, g=0, b=0; int k = 0;
106     while(k < 4 && u > knots[k+1]) k++;
107     m = (u-knots[k])/(knots[k+1]-knots[k]);
108     switch(k) {
109     case 0: r = 0;     g = 0;   b = (m+1)/2; break;
110     case 1: r = 0;     g = m;   b = 1;       break;
111     case 2: r = m;     g = 1;   b = 1-m;     break;
112     case 3: r = 1;     g = 1-m; b = 0;       break;
113     case 4: r = 1-m/2; g = 0;   b = 0;       break;
114     }
115     R[i] = (unsigned char)(255*PetscMin(r,1.0));
116     G[i] = (unsigned char)(255*PetscMin(g,1.0));
117     B[i] = (unsigned char)(255*PetscMin(b,1.0));
118   }
119   PetscFunctionReturn(0);
120 }
121 
122 #undef __FUNCT__
123 #define __FUNCT__ "PetscDrawCmap_Hot"
124 static PetscErrorCode PetscDrawCmap_Hot(int mapsize,unsigned char R[],unsigned char G[],unsigned char B[])
125 {
126   int          i;
127   const double knots[] =  {0, 3/8., 3/4., 1};
128 
129   PetscFunctionBegin;
130   for (i=0; i<mapsize; i++) {
131     double u = (double)i/(mapsize-1);
132     double m, r=0, g=0, b=0; int k = 0;
133     while(k < 2 && u > knots[k+1]) k++;
134     m = (u-knots[k])/(knots[k+1]-knots[k]);
135     switch(k) {
136     case 0: r = m; g = 0; b = 0; break;
137     case 1: r = 1; g = m; b = 0; break;
138     case 2: r = 1; g = 1; b = m; break;
139     }
140     R[i] = (unsigned char)(255*PetscMin(r,1.0));
141     G[i] = (unsigned char)(255*PetscMin(g,1.0));
142     B[i] = (unsigned char)(255*PetscMin(b,1.0));
143   }
144   PetscFunctionReturn(0);
145 }
146 
147 #undef __FUNCT__
148 #define __FUNCT__ "PetscDrawCmap_Bone"
149 static PetscErrorCode PetscDrawCmap_Bone(int mapsize,unsigned char R[],unsigned char G[],unsigned char B[])
150 {
151   int i;
152   PetscFunctionBegin;
153   (void)PetscDrawCmap_Hot(mapsize,R,G,B);
154   for (i=0; i<mapsize; i++) {
155     double u = (double)i/(mapsize-1);
156     double r = (7*u + B[i]/255.0)/8;
157     double g = (7*u + G[i]/255.0)/8;
158     double b = (7*u + R[i]/255.0)/8;
159     R[i] = (unsigned char)(255*PetscMin(r,1.0));
160     G[i] = (unsigned char)(255*PetscMin(g,1.0));
161     B[i] = (unsigned char)(255*PetscMin(b,1.0));
162   }
163   PetscFunctionReturn(0);
164 }
165 
166 #include "cmap/coolwarm.h"
167 #include "cmap/parula.h"
168 #include "cmap/viridis.h"
169 #include "cmap/plasma.h"
170 #include "cmap/inferno.h"
171 #include "cmap/magma.h"
172 
173 static struct {
174   const char           *name;
175   const unsigned char (*data)[3];
176   PetscErrorCode      (*cmap)(int,unsigned char[],unsigned char[],unsigned char[]);
177 } PetscDrawCmapTable[] = {
178   {"hue",      NULL, PetscDrawCmap_Hue },     /* varying hue with constant lightness and saturation */
179   {"gray",     NULL, PetscDrawCmap_Gray},     /* black to white with shades of gray */
180   {"bone",     NULL, PetscDrawCmap_Bone},     /* black to white with gray-blue shades */
181   {"jet",      NULL, PetscDrawCmap_Jet },     /* rainbow-like colormap from NCSA, University of Illinois */
182   {"hot",      NULL, PetscDrawCmap_Hot },     /* black-body radiation */
183   {"coolwarm", PetscDrawCmap_coolwarm, NULL}, /* ParaView default (Cool To Warm with Diverging interpolation) */
184   {"parula",   PetscDrawCmap_parula,   NULL}, /* MATLAB (default since R2014b) */
185   {"viridis",  PetscDrawCmap_viridis,  NULL}, /* matplotlib 1.5 (default since 2.0) */
186   {"plasma",   PetscDrawCmap_plasma,   NULL}, /* matplotlib 1.5 */
187   {"inferno",  PetscDrawCmap_inferno,  NULL}, /* matplotlib 1.5 */
188   {"magma",    PetscDrawCmap_magma,    NULL}, /* matplotlib 1.5 */
189 };
190 
191 #undef __FUNCT__
192 #define __FUNCT__ "PetscDrawUtilitySetCmap"
193 PetscErrorCode  PetscDrawUtilitySetCmap(const char colormap[],int mapsize,unsigned char R[],unsigned char G[],unsigned char B[])
194 {
195   int             i,j;
196   const char      *cmap_name_list[sizeof(PetscDrawCmapTable)/sizeof(PetscDrawCmapTable[0])];
197   PetscInt        id = 0, count = (PetscInt)(sizeof(cmap_name_list)/sizeof(char*));
198   PetscBool       reverse = PETSC_FALSE, brighten = PETSC_FALSE;
199   PetscReal       beta = 0;
200   PetscErrorCode  ierr;
201 
202   PetscFunctionBegin;
203   for (i=0; i<count; i++) cmap_name_list[i] = PetscDrawCmapTable[i].name;
204   if (colormap && colormap[0]) {
205     PetscBool match = PETSC_FALSE;
206     for (id=0; !match && id<count; id++) {ierr = PetscStrcasecmp(colormap,cmap_name_list[id],&match);CHKERRQ(ierr);}
207     if (!match) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Colormap '%s' not found",colormap);
208   }
209   ierr = PetscOptionsGetEList(NULL,NULL,"-draw_cmap",cmap_name_list,count,&id,NULL);CHKERRQ(ierr);
210   ierr = PetscOptionsGetBool(NULL,NULL,"-draw_cmap_reverse",&reverse,NULL);CHKERRQ(ierr);
211   ierr = PetscOptionsGetReal(NULL,NULL,"-draw_cmap_brighten",&beta,&brighten);CHKERRQ(ierr);
212   if (brighten && (beta <= (PetscReal)-1 || beta >= (PetscReal)+1)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"brighten parameter %g must be in the range (-1,1)",(double)beta);
213 
214   if (PetscDrawCmapTable[id].cmap) {
215     ierr = PetscDrawCmapTable[id].cmap(mapsize,R,G,B);CHKERRQ(ierr);
216   } else {
217     const unsigned char (*rgb)[3] = PetscDrawCmapTable[id].data;
218     if (mapsize != 256-PETSC_DRAW_BASIC_COLORS) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Colormap '%s' with size %d not supported",cmap_name_list[id],mapsize);
219     for (i=0; i<mapsize; i++) {R[i] = rgb[i][0]; G[i] = rgb[i][1]; B[i] = rgb[i][2];}
220   }
221 
222   if (reverse) {
223     i = 0; j = mapsize-1;
224     while(i < j) {
225 #define SWAP(a,i,j) do { unsigned char t = a[i]; a[i] = a[j]; a[j] = t; } while (0)
226       SWAP(R,i,j);
227       SWAP(G,i,j);
228       SWAP(B,i,j);
229 #undef SWAP
230       i++; j--;
231     }
232   }
233 
234   if (brighten) {
235     PetscReal gamma = (beta > 0.0) ? (1 - beta) : (1 / (1 + beta));
236     for (i=0; i<mapsize; i++) {
237       PetscReal r = PetscPowReal((PetscReal)R[i]/255,gamma);
238       PetscReal g = PetscPowReal((PetscReal)G[i]/255,gamma);
239       PetscReal b = PetscPowReal((PetscReal)B[i]/255,gamma);
240       R[i] = (unsigned char)(255*PetscMin(r,(PetscReal)1.0));
241       G[i] = (unsigned char)(255*PetscMin(g,(PetscReal)1.0));
242       B[i] = (unsigned char)(255*PetscMin(b,(PetscReal)1.0));
243     }
244   }
245   PetscFunctionReturn(0);
246 }
247