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 16d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDrawUtilitySetGamma(PetscReal g) 17d71ae5a4SJacob Faibussowitsch { 1800d931feSLisandro Dalcin PetscFunctionBegin; 1900d931feSLisandro Dalcin Gamma = g; 203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2100d931feSLisandro Dalcin } 2200d931feSLisandro Dalcin 23d71ae5a4SJacob Faibussowitsch static inline double PetscHlsHelper(double m1, double m2, double h) 24d71ae5a4SJacob Faibussowitsch { 2583b97956SLisandro Dalcin while (h > 1.0) h -= 1.0; 2683b97956SLisandro Dalcin while (h < 0.0) h += 1.0; 2783b97956SLisandro Dalcin if (h < 1 / 6.0) return m1 + (m2 - m1) * h * 6; 2883b97956SLisandro Dalcin if (h < 1 / 2.0) return m2; 2983b97956SLisandro Dalcin if (h < 2 / 3.0) return m1 + (m2 - m1) * (2 / 3.0 - h) * 6; 3000d931feSLisandro Dalcin return m1; 3100d931feSLisandro Dalcin } 3200d931feSLisandro Dalcin 33d71ae5a4SJacob Faibussowitsch static inline void PetscHlsToRgb(double h, double l, double s, double *r, double *g, double *b) 34d71ae5a4SJacob Faibussowitsch { 3583b97956SLisandro Dalcin if (s > 0.0) { 3683b97956SLisandro Dalcin double m2 = l <= 0.5 ? l * (1.0 + s) : l + s - (l * s); 3783b97956SLisandro Dalcin double m1 = 2 * l - m2; 3883b97956SLisandro Dalcin *r = PetscHlsHelper(m1, m2, h + 1 / 3.); 3983b97956SLisandro Dalcin *g = PetscHlsHelper(m1, m2, h); 4083b97956SLisandro Dalcin *b = PetscHlsHelper(m1, m2, h - 1 / 3.); 4100d931feSLisandro Dalcin } else { 4283b97956SLisandro Dalcin /* ignore hue */ 4383b97956SLisandro Dalcin *r = *g = *b = l; 4400d931feSLisandro Dalcin } 4500d931feSLisandro Dalcin } 4600d931feSLisandro Dalcin 47d71ae5a4SJacob Faibussowitsch static inline void PetscGammaCorrect(double *r, double *g, double *b) 48d71ae5a4SJacob Faibussowitsch { 4983b97956SLisandro Dalcin PetscReal igamma = 1 / Gamma; 5083b97956SLisandro Dalcin *r = (double)PetscPowReal((PetscReal)*r, igamma); 5183b97956SLisandro Dalcin *g = (double)PetscPowReal((PetscReal)*g, igamma); 5283b97956SLisandro Dalcin *b = (double)PetscPowReal((PetscReal)*b, igamma); 5383b97956SLisandro Dalcin } 5483b97956SLisandro Dalcin 55d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDrawCmap_Hue(int mapsize, unsigned char R[], unsigned char G[], unsigned char B[]) 56d71ae5a4SJacob Faibussowitsch { 5783b97956SLisandro Dalcin int i; 5883b97956SLisandro Dalcin double maxhue = 212.0 / 360, lightness = 0.5, saturation = 1.0; 5900d931feSLisandro Dalcin 6000d931feSLisandro Dalcin PetscFunctionBegin; 6100d931feSLisandro Dalcin for (i = 0; i < mapsize; i++) { 6283b97956SLisandro Dalcin double hue = maxhue * (double)i / (mapsize - 1), r, g, b; 6383b97956SLisandro Dalcin PetscHlsToRgb(hue, lightness, saturation, &r, &g, &b); 6483b97956SLisandro Dalcin PetscGammaCorrect(&r, &g, &b); 6583b97956SLisandro Dalcin R[i] = (unsigned char)(255 * PetscMin(r, 1.0)); 6683b97956SLisandro Dalcin G[i] = (unsigned char)(255 * PetscMin(g, 1.0)); 6783b97956SLisandro Dalcin B[i] = (unsigned char)(255 * PetscMin(b, 1.0)); 6800d931feSLisandro Dalcin } 693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7000d931feSLisandro Dalcin } 7100d931feSLisandro Dalcin 72d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDrawCmap_Gray(int mapsize, unsigned char R[], unsigned char G[], unsigned char B[]) 73d71ae5a4SJacob Faibussowitsch { 7400d931feSLisandro Dalcin PetscFunctionBegin; 75*4d86920dSPierre Jolivet for (int i = 0; i < mapsize; i++) R[i] = G[i] = B[i] = (unsigned char)((255.0 * i) / (mapsize - 1)); 763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7700d931feSLisandro Dalcin } 7800d931feSLisandro Dalcin 79d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDrawCmap_Jet(int mapsize, unsigned char R[], unsigned char G[], unsigned char B[]) 80d71ae5a4SJacob Faibussowitsch { 8100d931feSLisandro Dalcin int i; 8200d931feSLisandro Dalcin const double knots[] = {0, 1 / 8., 3 / 8., 5 / 8., 7 / 8., 1}; 8300d931feSLisandro Dalcin 8400d931feSLisandro Dalcin PetscFunctionBegin; 8500d931feSLisandro Dalcin for (i = 0; i < mapsize; i++) { 8600d931feSLisandro Dalcin double u = (double)i / (mapsize - 1); 879371c9d4SSatish Balay double m, r = 0, g = 0, b = 0; 889371c9d4SSatish Balay int k = 0; 8900d931feSLisandro Dalcin while (k < 4 && u > knots[k + 1]) k++; 9000d931feSLisandro Dalcin m = (u - knots[k]) / (knots[k + 1] - knots[k]); 9100d931feSLisandro Dalcin switch (k) { 929371c9d4SSatish Balay case 0: 939371c9d4SSatish Balay r = 0; 949371c9d4SSatish Balay g = 0; 959371c9d4SSatish Balay b = (m + 1) / 2; 969371c9d4SSatish Balay break; 979371c9d4SSatish Balay case 1: 989371c9d4SSatish Balay r = 0; 999371c9d4SSatish Balay g = m; 1009371c9d4SSatish Balay b = 1; 1019371c9d4SSatish Balay break; 1029371c9d4SSatish Balay case 2: 1039371c9d4SSatish Balay r = m; 1049371c9d4SSatish Balay g = 1; 1059371c9d4SSatish Balay b = 1 - m; 1069371c9d4SSatish Balay break; 1079371c9d4SSatish Balay case 3: 1089371c9d4SSatish Balay r = 1; 1099371c9d4SSatish Balay g = 1 - m; 1109371c9d4SSatish Balay b = 0; 1119371c9d4SSatish Balay break; 1129371c9d4SSatish Balay case 4: 1139371c9d4SSatish Balay r = 1 - m / 2; 1149371c9d4SSatish Balay g = 0; 1159371c9d4SSatish Balay b = 0; 1169371c9d4SSatish Balay break; 11700d931feSLisandro Dalcin } 11800d931feSLisandro Dalcin R[i] = (unsigned char)(255 * PetscMin(r, 1.0)); 11900d931feSLisandro Dalcin G[i] = (unsigned char)(255 * PetscMin(g, 1.0)); 12000d931feSLisandro Dalcin B[i] = (unsigned char)(255 * PetscMin(b, 1.0)); 12100d931feSLisandro Dalcin } 1223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12300d931feSLisandro Dalcin } 12400d931feSLisandro Dalcin 125d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDrawCmap_Hot(int mapsize, unsigned char R[], unsigned char G[], unsigned char B[]) 126d71ae5a4SJacob Faibussowitsch { 12700d931feSLisandro Dalcin int i; 12800d931feSLisandro Dalcin const double knots[] = {0, 3 / 8., 3 / 4., 1}; 12900d931feSLisandro Dalcin 13000d931feSLisandro Dalcin PetscFunctionBegin; 13100d931feSLisandro Dalcin for (i = 0; i < mapsize; i++) { 13200d931feSLisandro Dalcin double u = (double)i / (mapsize - 1); 1339371c9d4SSatish Balay double m, r = 0, g = 0, b = 0; 1349371c9d4SSatish Balay int k = 0; 13500d931feSLisandro Dalcin while (k < 2 && u > knots[k + 1]) k++; 13600d931feSLisandro Dalcin m = (u - knots[k]) / (knots[k + 1] - knots[k]); 13700d931feSLisandro Dalcin switch (k) { 1389371c9d4SSatish Balay case 0: 1399371c9d4SSatish Balay r = m; 1409371c9d4SSatish Balay g = 0; 1419371c9d4SSatish Balay b = 0; 1429371c9d4SSatish Balay break; 1439371c9d4SSatish Balay case 1: 1449371c9d4SSatish Balay r = 1; 1459371c9d4SSatish Balay g = m; 1469371c9d4SSatish Balay b = 0; 1479371c9d4SSatish Balay break; 1489371c9d4SSatish Balay case 2: 1499371c9d4SSatish Balay r = 1; 1509371c9d4SSatish Balay g = 1; 1519371c9d4SSatish Balay b = m; 1529371c9d4SSatish Balay break; 15300d931feSLisandro Dalcin } 15400d931feSLisandro Dalcin R[i] = (unsigned char)(255 * PetscMin(r, 1.0)); 15500d931feSLisandro Dalcin G[i] = (unsigned char)(255 * PetscMin(g, 1.0)); 15600d931feSLisandro Dalcin B[i] = (unsigned char)(255 * PetscMin(b, 1.0)); 15700d931feSLisandro Dalcin } 1583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15900d931feSLisandro Dalcin } 16000d931feSLisandro Dalcin 161d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDrawCmap_Bone(int mapsize, unsigned char R[], unsigned char G[], unsigned char B[]) 162d71ae5a4SJacob Faibussowitsch { 16300d931feSLisandro Dalcin PetscFunctionBegin; 16400d931feSLisandro Dalcin (void)PetscDrawCmap_Hot(mapsize, R, G, B); 165*4d86920dSPierre Jolivet for (int i = 0; i < mapsize; i++) { 16600d931feSLisandro Dalcin double u = (double)i / (mapsize - 1); 16700d931feSLisandro Dalcin double r = (7 * u + B[i] / 255.0) / 8; 16800d931feSLisandro Dalcin double g = (7 * u + G[i] / 255.0) / 8; 16900d931feSLisandro Dalcin double b = (7 * u + R[i] / 255.0) / 8; 17000d931feSLisandro Dalcin R[i] = (unsigned char)(255 * PetscMin(r, 1.0)); 17100d931feSLisandro Dalcin G[i] = (unsigned char)(255 * PetscMin(g, 1.0)); 17200d931feSLisandro Dalcin B[i] = (unsigned char)(255 * PetscMin(b, 1.0)); 17300d931feSLisandro Dalcin } 1743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17500d931feSLisandro Dalcin } 17600d931feSLisandro Dalcin 177575a0592SBarry Smith #include "../src/sys/classes/draw/utils/cmap/coolwarm.h" 178575a0592SBarry Smith #include "../src/sys/classes/draw/utils/cmap/parula.h" 179575a0592SBarry Smith #include "../src/sys/classes/draw/utils/cmap/viridis.h" 180575a0592SBarry Smith #include "../src/sys/classes/draw/utils/cmap/plasma.h" 181575a0592SBarry Smith #include "../src/sys/classes/draw/utils/cmap/inferno.h" 182575a0592SBarry Smith #include "../src/sys/classes/draw/utils/cmap/magma.h" 18300d931feSLisandro Dalcin 18400d931feSLisandro Dalcin static struct { 18500d931feSLisandro Dalcin const char *name; 18600d931feSLisandro Dalcin const unsigned char (*data)[3]; 18700d931feSLisandro Dalcin PetscErrorCode (*cmap)(int, unsigned char[], unsigned char[], unsigned char[]); 18800d931feSLisandro Dalcin } PetscDrawCmapTable[] = { 18900d931feSLisandro Dalcin {"hue", NULL, PetscDrawCmap_Hue }, /* varying hue with constant lightness and saturation */ 19000d931feSLisandro Dalcin {"gray", NULL, PetscDrawCmap_Gray}, /* black to white with shades of gray */ 19100d931feSLisandro Dalcin {"bone", NULL, PetscDrawCmap_Bone}, /* black to white with gray-blue shades */ 19200d931feSLisandro Dalcin {"jet", NULL, PetscDrawCmap_Jet }, /* rainbow-like colormap from NCSA, University of Illinois */ 19300d931feSLisandro Dalcin {"hot", NULL, PetscDrawCmap_Hot }, /* black-body radiation */ 19400d931feSLisandro Dalcin {"coolwarm", PetscDrawCmap_coolwarm, NULL }, /* ParaView default (Cool To Warm with Diverging interpolation) */ 19500d931feSLisandro Dalcin {"parula", PetscDrawCmap_parula, NULL }, /* MATLAB (default since R2014b) */ 19600d931feSLisandro Dalcin {"viridis", PetscDrawCmap_viridis, NULL }, /* matplotlib 1.5 (default since 2.0) */ 19700d931feSLisandro Dalcin {"plasma", PetscDrawCmap_plasma, NULL }, /* matplotlib 1.5 */ 19800d931feSLisandro Dalcin {"inferno", PetscDrawCmap_inferno, NULL }, /* matplotlib 1.5 */ 19900d931feSLisandro Dalcin {"magma", PetscDrawCmap_magma, NULL }, /* matplotlib 1.5 */ 20000d931feSLisandro Dalcin }; 20100d931feSLisandro Dalcin 202d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDrawUtilitySetCmap(const char colormap[], int mapsize, unsigned char R[], unsigned char G[], unsigned char B[]) 203d71ae5a4SJacob Faibussowitsch { 20400d931feSLisandro Dalcin int i, j; 205dd39110bSPierre Jolivet const char *cmap_name_list[PETSC_STATIC_ARRAY_LENGTH(PetscDrawCmapTable)]; 206400f6f02SBarry Smith PetscInt id = 0, count = PETSC_STATIC_ARRAY_LENGTH(cmap_name_list); 20700d931feSLisandro Dalcin PetscBool reverse = PETSC_FALSE, brighten = PETSC_FALSE; 20800d931feSLisandro Dalcin PetscReal beta = 0; 20900d931feSLisandro Dalcin 21000d931feSLisandro Dalcin PetscFunctionBegin; 21100d931feSLisandro Dalcin for (i = 0; i < count; i++) cmap_name_list[i] = PetscDrawCmapTable[i].name; 21200d931feSLisandro Dalcin if (colormap && colormap[0]) { 21300d931feSLisandro Dalcin PetscBool match = PETSC_FALSE; 2149566063dSJacob Faibussowitsch for (id = 0; !match && id < count; id++) PetscCall(PetscStrcasecmp(colormap, cmap_name_list[id], &match)); 21528b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Colormap '%s' not found", colormap); 21600d931feSLisandro Dalcin } 2179566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEList(NULL, NULL, "-draw_cmap", cmap_name_list, count, &id, NULL)); 2189566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_cmap_reverse", &reverse, NULL)); 2199566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-draw_cmap_brighten", &beta, &brighten)); 220cc73adaaSBarry Smith PetscCheck(!brighten || (beta > (PetscReal)-1 && beta < (PetscReal) + 1), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "brighten parameter %g must be in the range (-1,1)", (double)beta); 22100d931feSLisandro Dalcin 22200d931feSLisandro Dalcin if (PetscDrawCmapTable[id].cmap) { 2239566063dSJacob Faibussowitsch PetscCall(PetscDrawCmapTable[id].cmap(mapsize, R, G, B)); 22400d931feSLisandro Dalcin } else { 22500d931feSLisandro Dalcin const unsigned char(*rgb)[3] = PetscDrawCmapTable[id].data; 22608401ef6SPierre Jolivet PetscCheck(mapsize == 256 - PETSC_DRAW_BASIC_COLORS, PETSC_COMM_SELF, PETSC_ERR_SUP, "Colormap '%s' with size %d not supported", cmap_name_list[id], mapsize); 2279371c9d4SSatish Balay for (i = 0; i < mapsize; i++) { 2289371c9d4SSatish Balay R[i] = rgb[i][0]; 2299371c9d4SSatish Balay G[i] = rgb[i][1]; 2309371c9d4SSatish Balay B[i] = rgb[i][2]; 2319371c9d4SSatish Balay } 23200d931feSLisandro Dalcin } 23300d931feSLisandro Dalcin 23400d931feSLisandro Dalcin if (reverse) { 2359371c9d4SSatish Balay i = 0; 2369371c9d4SSatish Balay j = mapsize - 1; 23700d931feSLisandro Dalcin while (i < j) { 2389371c9d4SSatish Balay #define SWAP(a, i, j) \ 2399371c9d4SSatish Balay do { \ 2409371c9d4SSatish Balay unsigned char t = a[i]; \ 2419371c9d4SSatish Balay a[i] = a[j]; \ 2429371c9d4SSatish Balay a[j] = t; \ 2439371c9d4SSatish Balay } while (0) 24400d931feSLisandro Dalcin SWAP(R, i, j); 24500d931feSLisandro Dalcin SWAP(G, i, j); 24600d931feSLisandro Dalcin SWAP(B, i, j); 24700d931feSLisandro Dalcin #undef SWAP 2489371c9d4SSatish Balay i++; 2499371c9d4SSatish Balay j--; 25000d931feSLisandro Dalcin } 25100d931feSLisandro Dalcin } 25200d931feSLisandro Dalcin 25300d931feSLisandro Dalcin if (brighten) { 25400d931feSLisandro Dalcin PetscReal gamma = (beta > 0.0) ? (1 - beta) : (1 / (1 + beta)); 25500d931feSLisandro Dalcin for (i = 0; i < mapsize; i++) { 25600d931feSLisandro Dalcin PetscReal r = PetscPowReal((PetscReal)R[i] / 255, gamma); 25700d931feSLisandro Dalcin PetscReal g = PetscPowReal((PetscReal)G[i] / 255, gamma); 25800d931feSLisandro Dalcin PetscReal b = PetscPowReal((PetscReal)B[i] / 255, gamma); 25900d931feSLisandro Dalcin R[i] = (unsigned char)(255 * PetscMin(r, (PetscReal)1.0)); 26000d931feSLisandro Dalcin G[i] = (unsigned char)(255 * PetscMin(g, (PetscReal)1.0)); 26100d931feSLisandro Dalcin B[i] = (unsigned char)(255 * PetscMin(b, (PetscReal)1.0)); 26200d931feSLisandro Dalcin } 26300d931feSLisandro Dalcin } 2643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 26500d931feSLisandro Dalcin } 266