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