xref: /petsc/src/sys/classes/draw/utils/cmap.c (revision 83b979563e1b6cd021d59e9811fc6fa38b93af5f)
100d931feSLisandro Dalcin #include <petscsys.h>              /*I "petscsys.h" I*/
200d931feSLisandro Dalcin #include <petscdraw.h>
300d931feSLisandro Dalcin 
400d931feSLisandro Dalcin /*
500d931feSLisandro Dalcin     Set up a color map, using uniform separation in hue space.
600d931feSLisandro Dalcin     Map entries are Red, Green, Blue.
700d931feSLisandro Dalcin     Values are "gamma" corrected.
800d931feSLisandro Dalcin  */
900d931feSLisandro Dalcin 
1000d931feSLisandro Dalcin /*
1100d931feSLisandro Dalcin    Gamma is a monitor dependent value.  The value here is an
1200d931feSLisandro Dalcin    approximate that gives somewhat better results than Gamma = 1.
1300d931feSLisandro Dalcin  */
1400d931feSLisandro Dalcin static PetscReal Gamma = 2.0;
1500d931feSLisandro Dalcin 
1600d931feSLisandro Dalcin #undef __FUNCT__
1700d931feSLisandro Dalcin #define __FUNCT__ "PetscDrawUtilitySetGamma"
1800d931feSLisandro Dalcin PetscErrorCode  PetscDrawUtilitySetGamma(PetscReal g)
1900d931feSLisandro Dalcin {
2000d931feSLisandro Dalcin   PetscFunctionBegin;
2100d931feSLisandro Dalcin   Gamma = g;
2200d931feSLisandro Dalcin   PetscFunctionReturn(0);
2300d931feSLisandro Dalcin }
2400d931feSLisandro Dalcin 
25*83b97956SLisandro Dalcin PETSC_STATIC_INLINE double PetscHlsHelper(double m1,double m2,double h)
2600d931feSLisandro Dalcin {
27*83b97956SLisandro Dalcin   while (h > 1.0) h -= 1.0;
28*83b97956SLisandro Dalcin   while (h < 0.0) h += 1.0;
29*83b97956SLisandro Dalcin   if (h < 1/6.0) return m1 + (m2-m1)*h*6;
30*83b97956SLisandro Dalcin   if (h < 1/2.0) return m2;
31*83b97956SLisandro Dalcin   if (h < 2/3.0) return m1 + (m2-m1)*(2/3.0-h)*6;
3200d931feSLisandro Dalcin   return m1;
3300d931feSLisandro Dalcin }
3400d931feSLisandro Dalcin 
35*83b97956SLisandro Dalcin PETSC_STATIC_INLINE void PetscHlsToRgb(double h,double l,double s,double *r,double *g,double *b)
3600d931feSLisandro Dalcin {
37*83b97956SLisandro Dalcin   if (s > 0.0) {
38*83b97956SLisandro Dalcin     double m2 = l <= 0.5 ? l * (1.0+s) : l+s-(l*s);
39*83b97956SLisandro Dalcin     double m1 = 2*l - m2;
40*83b97956SLisandro Dalcin     *r = PetscHlsHelper(m1,m2,h+1/3.);
41*83b97956SLisandro Dalcin     *g = PetscHlsHelper(m1,m2,h);
42*83b97956SLisandro Dalcin     *b = PetscHlsHelper(m1,m2,h-1/3.);
4300d931feSLisandro Dalcin   } else {
44*83b97956SLisandro Dalcin     /* ignore hue */
45*83b97956SLisandro Dalcin     *r = *g = *b = l;
4600d931feSLisandro Dalcin   }
4700d931feSLisandro Dalcin }
4800d931feSLisandro Dalcin 
49*83b97956SLisandro Dalcin PETSC_STATIC_INLINE void PetscGammaCorrect(double *r,double *g,double *b)
50*83b97956SLisandro Dalcin {
51*83b97956SLisandro Dalcin   PetscReal igamma = 1/Gamma;
52*83b97956SLisandro Dalcin   *r = (double)PetscPowReal((PetscReal)*r,igamma);
53*83b97956SLisandro Dalcin   *g = (double)PetscPowReal((PetscReal)*g,igamma);
54*83b97956SLisandro Dalcin   *b = (double)PetscPowReal((PetscReal)*b,igamma);
55*83b97956SLisandro Dalcin }
56*83b97956SLisandro Dalcin 
5700d931feSLisandro Dalcin #undef __FUNCT__
5800d931feSLisandro Dalcin #define __FUNCT__ "PetscDrawCmap_Hue"
5900d931feSLisandro Dalcin static PetscErrorCode PetscDrawCmap_Hue(int mapsize, unsigned char R[],unsigned char G[],unsigned char B[])
6000d931feSLisandro Dalcin {
61*83b97956SLisandro Dalcin   int    i;
62*83b97956SLisandro Dalcin   double maxhue = 212.0/360,lightness = 0.5,saturation = 1.0;
6300d931feSLisandro Dalcin 
6400d931feSLisandro Dalcin   PetscFunctionBegin;
6500d931feSLisandro Dalcin   for (i=0; i<mapsize; i++) {
66*83b97956SLisandro Dalcin     double hue = maxhue*(double)i/(mapsize-1),r,g,b;
67*83b97956SLisandro Dalcin     PetscHlsToRgb(hue,lightness,saturation,&r,&g,&b);
68*83b97956SLisandro Dalcin     PetscGammaCorrect(&r,&g,&b);
69*83b97956SLisandro Dalcin     R[i] = (unsigned char)(255*PetscMin(r,1.0));
70*83b97956SLisandro Dalcin     G[i] = (unsigned char)(255*PetscMin(g,1.0));
71*83b97956SLisandro Dalcin     B[i] = (unsigned char)(255*PetscMin(b,1.0));
7200d931feSLisandro Dalcin   }
7300d931feSLisandro Dalcin   PetscFunctionReturn(0);
7400d931feSLisandro Dalcin }
7500d931feSLisandro Dalcin 
7600d931feSLisandro Dalcin #undef __FUNCT__
7700d931feSLisandro Dalcin #define __FUNCT__ "PetscDrawCmap_Gray"
7800d931feSLisandro Dalcin static PetscErrorCode PetscDrawCmap_Gray(int mapsize,unsigned char R[],unsigned char G[],unsigned char B[])
7900d931feSLisandro Dalcin {
8000d931feSLisandro Dalcin   int i;
8100d931feSLisandro Dalcin   PetscFunctionBegin;
8200d931feSLisandro Dalcin   for (i=0; i<mapsize; i++) R[i] = G[i] = B[i] = (unsigned char)((255.0*i)/(mapsize-1));
8300d931feSLisandro Dalcin   PetscFunctionReturn(0);
8400d931feSLisandro Dalcin }
8500d931feSLisandro Dalcin 
8600d931feSLisandro Dalcin #undef __FUNCT__
8700d931feSLisandro Dalcin #define __FUNCT__ "PetscDrawCmap_Jet"
8800d931feSLisandro Dalcin static PetscErrorCode PetscDrawCmap_Jet(int mapsize,unsigned char R[],unsigned char G[],unsigned char B[])
8900d931feSLisandro Dalcin {
9000d931feSLisandro Dalcin   int          i;
9100d931feSLisandro Dalcin   const double knots[] =  {0, 1/8., 3/8., 5/8., 7/8., 1};
9200d931feSLisandro Dalcin 
9300d931feSLisandro Dalcin   PetscFunctionBegin;
9400d931feSLisandro Dalcin   for (i=0; i<mapsize; i++) {
9500d931feSLisandro Dalcin     double u = (double)i/(mapsize-1);
9600d931feSLisandro Dalcin     double m, r=0, g=0, b=0; int k = 0;
9700d931feSLisandro Dalcin     while(k < 4 && u > knots[k+1]) k++;
9800d931feSLisandro Dalcin     m = (u-knots[k])/(knots[k+1]-knots[k]);
9900d931feSLisandro Dalcin     switch(k) {
10000d931feSLisandro Dalcin     case 0: r = 0;     g = 0;   b = (m+1)/2; break;
10100d931feSLisandro Dalcin     case 1: r = 0;     g = m;   b = 1;       break;
10200d931feSLisandro Dalcin     case 2: r = m;     g = 1;   b = 1-m;     break;
10300d931feSLisandro Dalcin     case 3: r = 1;     g = 1-m; b = 0;       break;
10400d931feSLisandro Dalcin     case 4: r = 1-m/2; g = 0;   b = 0;       break;
10500d931feSLisandro Dalcin     }
10600d931feSLisandro Dalcin     R[i] = (unsigned char)(255*PetscMin(r,1.0));
10700d931feSLisandro Dalcin     G[i] = (unsigned char)(255*PetscMin(g,1.0));
10800d931feSLisandro Dalcin     B[i] = (unsigned char)(255*PetscMin(b,1.0));
10900d931feSLisandro Dalcin   }
11000d931feSLisandro Dalcin   PetscFunctionReturn(0);
11100d931feSLisandro Dalcin }
11200d931feSLisandro Dalcin 
11300d931feSLisandro Dalcin #undef __FUNCT__
11400d931feSLisandro Dalcin #define __FUNCT__ "PetscDrawCmap_Hot"
11500d931feSLisandro Dalcin static PetscErrorCode PetscDrawCmap_Hot(int mapsize,unsigned char R[],unsigned char G[],unsigned char B[])
11600d931feSLisandro Dalcin {
11700d931feSLisandro Dalcin   int          i;
11800d931feSLisandro Dalcin   const double knots[] =  {0, 3/8., 3/4., 1};
11900d931feSLisandro Dalcin 
12000d931feSLisandro Dalcin   PetscFunctionBegin;
12100d931feSLisandro Dalcin   for (i=0; i<mapsize; i++) {
12200d931feSLisandro Dalcin     double u = (double)i/(mapsize-1);
12300d931feSLisandro Dalcin     double m, r=0, g=0, b=0; int k = 0;
12400d931feSLisandro Dalcin     while(k < 2 && u > knots[k+1]) k++;
12500d931feSLisandro Dalcin     m = (u-knots[k])/(knots[k+1]-knots[k]);
12600d931feSLisandro Dalcin     switch(k) {
12700d931feSLisandro Dalcin     case 0: r = m; g = 0; b = 0; break;
12800d931feSLisandro Dalcin     case 1: r = 1; g = m; b = 0; break;
12900d931feSLisandro Dalcin     case 2: r = 1; g = 1; b = m; break;
13000d931feSLisandro Dalcin     }
13100d931feSLisandro Dalcin     R[i] = (unsigned char)(255*PetscMin(r,1.0));
13200d931feSLisandro Dalcin     G[i] = (unsigned char)(255*PetscMin(g,1.0));
13300d931feSLisandro Dalcin     B[i] = (unsigned char)(255*PetscMin(b,1.0));
13400d931feSLisandro Dalcin   }
13500d931feSLisandro Dalcin   PetscFunctionReturn(0);
13600d931feSLisandro Dalcin }
13700d931feSLisandro Dalcin 
13800d931feSLisandro Dalcin #undef __FUNCT__
13900d931feSLisandro Dalcin #define __FUNCT__ "PetscDrawCmap_Bone"
14000d931feSLisandro Dalcin static PetscErrorCode PetscDrawCmap_Bone(int mapsize,unsigned char R[],unsigned char G[],unsigned char B[])
14100d931feSLisandro Dalcin {
14200d931feSLisandro Dalcin   int i;
14300d931feSLisandro Dalcin   PetscFunctionBegin;
14400d931feSLisandro Dalcin   (void)PetscDrawCmap_Hot(mapsize,R,G,B);
14500d931feSLisandro Dalcin   for (i=0; i<mapsize; i++) {
14600d931feSLisandro Dalcin     double u = (double)i/(mapsize-1);
14700d931feSLisandro Dalcin     double r = (7*u + B[i]/255.0)/8;
14800d931feSLisandro Dalcin     double g = (7*u + G[i]/255.0)/8;
14900d931feSLisandro Dalcin     double b = (7*u + R[i]/255.0)/8;
15000d931feSLisandro Dalcin     R[i] = (unsigned char)(255*PetscMin(r,1.0));
15100d931feSLisandro Dalcin     G[i] = (unsigned char)(255*PetscMin(g,1.0));
15200d931feSLisandro Dalcin     B[i] = (unsigned char)(255*PetscMin(b,1.0));
15300d931feSLisandro Dalcin   }
15400d931feSLisandro Dalcin   PetscFunctionReturn(0);
15500d931feSLisandro Dalcin }
15600d931feSLisandro Dalcin 
15700d931feSLisandro Dalcin #include "cmap/coolwarm.h"
15800d931feSLisandro Dalcin #include "cmap/parula.h"
15900d931feSLisandro Dalcin #include "cmap/viridis.h"
16000d931feSLisandro Dalcin #include "cmap/plasma.h"
16100d931feSLisandro Dalcin #include "cmap/inferno.h"
16200d931feSLisandro Dalcin #include "cmap/magma.h"
16300d931feSLisandro Dalcin 
16400d931feSLisandro Dalcin static struct {
16500d931feSLisandro Dalcin   const char           *name;
16600d931feSLisandro Dalcin   const unsigned char (*data)[3];
16700d931feSLisandro Dalcin   PetscErrorCode      (*cmap)(int,unsigned char[],unsigned char[],unsigned char[]);
16800d931feSLisandro Dalcin } PetscDrawCmapTable[] = {
16900d931feSLisandro Dalcin   {"hue",      NULL, PetscDrawCmap_Hue },     /* varying hue with constant lightness and saturation */
17000d931feSLisandro Dalcin   {"gray",     NULL, PetscDrawCmap_Gray},     /* black to white with shades of gray */
17100d931feSLisandro Dalcin   {"bone",     NULL, PetscDrawCmap_Bone},     /* black to white with gray-blue shades */
17200d931feSLisandro Dalcin   {"jet",      NULL, PetscDrawCmap_Jet },     /* rainbow-like colormap from NCSA, University of Illinois */
17300d931feSLisandro Dalcin   {"hot",      NULL, PetscDrawCmap_Hot },     /* black-body radiation */
17400d931feSLisandro Dalcin   {"coolwarm", PetscDrawCmap_coolwarm, NULL}, /* ParaView default (Cool To Warm with Diverging interpolation) */
17500d931feSLisandro Dalcin   {"parula",   PetscDrawCmap_parula,   NULL}, /* MATLAB (default since R2014b) */
17600d931feSLisandro Dalcin   {"viridis",  PetscDrawCmap_viridis,  NULL}, /* matplotlib 1.5 (default since 2.0) */
17700d931feSLisandro Dalcin   {"plasma",   PetscDrawCmap_plasma,   NULL}, /* matplotlib 1.5 */
17800d931feSLisandro Dalcin   {"inferno",  PetscDrawCmap_inferno,  NULL}, /* matplotlib 1.5 */
17900d931feSLisandro Dalcin   {"magma",    PetscDrawCmap_magma,    NULL}, /* matplotlib 1.5 */
18000d931feSLisandro Dalcin };
18100d931feSLisandro Dalcin 
18200d931feSLisandro Dalcin #undef __FUNCT__
18300d931feSLisandro Dalcin #define __FUNCT__ "PetscDrawUtilitySetCmap"
18400d931feSLisandro Dalcin PetscErrorCode  PetscDrawUtilitySetCmap(const char colormap[],int mapsize,unsigned char R[],unsigned char G[],unsigned char B[])
18500d931feSLisandro Dalcin {
18600d931feSLisandro Dalcin   int             i,j;
18700d931feSLisandro Dalcin   const char      *cmap_name_list[sizeof(PetscDrawCmapTable)/sizeof(PetscDrawCmapTable[0])];
18800d931feSLisandro Dalcin   PetscInt        id = 0, count = (PetscInt)(sizeof(cmap_name_list)/sizeof(char*));
18900d931feSLisandro Dalcin   PetscBool       reverse = PETSC_FALSE, brighten = PETSC_FALSE;
19000d931feSLisandro Dalcin   PetscReal       beta = 0;
19100d931feSLisandro Dalcin   PetscErrorCode  ierr;
19200d931feSLisandro Dalcin 
19300d931feSLisandro Dalcin   PetscFunctionBegin;
19400d931feSLisandro Dalcin   for (i=0; i<count; i++) cmap_name_list[i] = PetscDrawCmapTable[i].name;
19500d931feSLisandro Dalcin   if (colormap && colormap[0]) {
19600d931feSLisandro Dalcin     PetscBool match = PETSC_FALSE;
19700d931feSLisandro Dalcin     for (id=0; !match && id<count; id++) {ierr = PetscStrcasecmp(colormap,cmap_name_list[id],&match);CHKERRQ(ierr);}
19800d931feSLisandro Dalcin     if (!match) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Colormap '%s' not found",colormap);
19900d931feSLisandro Dalcin   }
20000d931feSLisandro Dalcin   ierr = PetscOptionsGetEList(NULL,NULL,"-draw_cmap",cmap_name_list,count,&id,NULL);CHKERRQ(ierr);
20100d931feSLisandro Dalcin   ierr = PetscOptionsGetBool(NULL,NULL,"-draw_cmap_reverse",&reverse,NULL);CHKERRQ(ierr);
20200d931feSLisandro Dalcin   ierr = PetscOptionsGetReal(NULL,NULL,"-draw_cmap_brighten",&beta,&brighten);CHKERRQ(ierr);
20300d931feSLisandro 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);
20400d931feSLisandro Dalcin 
20500d931feSLisandro Dalcin   if (PetscDrawCmapTable[id].cmap) {
20600d931feSLisandro Dalcin     ierr = PetscDrawCmapTable[id].cmap(mapsize,R,G,B);CHKERRQ(ierr);
20700d931feSLisandro Dalcin   } else {
20800d931feSLisandro Dalcin     const unsigned char (*rgb)[3] = PetscDrawCmapTable[id].data;
20900d931feSLisandro 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);
21000d931feSLisandro Dalcin     for (i=0; i<mapsize; i++) {R[i] = rgb[i][0]; G[i] = rgb[i][1]; B[i] = rgb[i][2];}
21100d931feSLisandro Dalcin   }
21200d931feSLisandro Dalcin 
21300d931feSLisandro Dalcin   if (reverse) {
21400d931feSLisandro Dalcin     i = 0; j = mapsize-1;
21500d931feSLisandro Dalcin     while(i < j) {
21600d931feSLisandro Dalcin #define SWAP(a,i,j) do { unsigned char t = a[i]; a[i] = a[j]; a[j] = t; } while (0)
21700d931feSLisandro Dalcin       SWAP(R,i,j);
21800d931feSLisandro Dalcin       SWAP(G,i,j);
21900d931feSLisandro Dalcin       SWAP(B,i,j);
22000d931feSLisandro Dalcin #undef SWAP
22100d931feSLisandro Dalcin       i++; j--;
22200d931feSLisandro Dalcin     }
22300d931feSLisandro Dalcin   }
22400d931feSLisandro Dalcin 
22500d931feSLisandro Dalcin   if (brighten) {
22600d931feSLisandro Dalcin     PetscReal gamma = (beta > 0.0) ? (1 - beta) : (1 / (1 + beta));
22700d931feSLisandro Dalcin     for (i=0; i<mapsize; i++) {
22800d931feSLisandro Dalcin       PetscReal r = PetscPowReal((PetscReal)R[i]/255,gamma);
22900d931feSLisandro Dalcin       PetscReal g = PetscPowReal((PetscReal)G[i]/255,gamma);
23000d931feSLisandro Dalcin       PetscReal b = PetscPowReal((PetscReal)B[i]/255,gamma);
23100d931feSLisandro Dalcin       R[i] = (unsigned char)(255*PetscMin(r,(PetscReal)1.0));
23200d931feSLisandro Dalcin       G[i] = (unsigned char)(255*PetscMin(g,(PetscReal)1.0));
23300d931feSLisandro Dalcin       B[i] = (unsigned char)(255*PetscMin(b,(PetscReal)1.0));
23400d931feSLisandro Dalcin     }
23500d931feSLisandro Dalcin   }
23600d931feSLisandro Dalcin   PetscFunctionReturn(0);
23700d931feSLisandro Dalcin }
238