Logo ROOT   master
Reference Guide
TColor.cxx
Go to the documentation of this file.
1 // @(#)root/base:$Id$
2 // Author: Rene Brun 12/12/94
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 #include "TROOT.h"
13 #include "TColor.h"
14 #include "TObjArray.h"
15 #include "TArrayI.h"
16 #include "TArrayD.h"
17 #include "TVirtualX.h"
18 #include "TError.h"
19 #include "TMathBase.h"
20 #include "TApplication.h"
21 #include <algorithm>
22 #include <cmath>
23 #include <iostream>
24 
26 
27 namespace {
28  static Bool_t& TColor__GrayScaleMode() {
29  static Bool_t grayScaleMode;
30  return grayScaleMode;
31  }
32  static TArrayI& TColor__Palette() {
33  static TArrayI globalPalette(0);
34  return globalPalette;
35  }
36  static TArrayD& TColor__PalettesList() {
37  static TArrayD globalPalettesList(0);
38  return globalPalettesList;
39  }
40 }
41 
42 static Int_t gHighestColorIndex = 0; ///< Highest color index defined
43 static Float_t gColorThreshold = -1.; ///< Color threshold used by GetColor
44 static Int_t gDefinedColors = 0; ///< Number of defined colors.
45 static Int_t gLastDefinedColors = 649; ///< Previous number of defined colors
46 
47 #define fgGrayscaleMode TColor__GrayScaleMode()
48 #define fgPalette TColor__Palette()
49 #define fgPalettesList TColor__PalettesList()
50 
51 using std::floor;
52 
53 /** \class TColor
54 \ingroup Base
55 \ingroup GraphicsAtt
56 
57 The color creation and management class.
58 
59  - [Introduction](#C00)
60  - [Basic colors](#C01)
61  - [The color wheel](#C02)
62  - [Bright and dark colors](#C03)
63  - [Gray scale view of of canvas with colors](#C04)
64  - [Color palettes](#C05)
65  - [High quality predefined palettes](#C06)
66  - [Palette inversion](#C061)
67  - [Color transparency](#C07)
68 
69 ## <a name="C00"></a> Introduction
70 
71 Colors are defined by their red, green and blue components, simply called the
72 RGB components. The colors are also known by the hue, light and saturation
73 components also known as the HLS components. When a new color is created the
74 components of both color systems are computed.
75 
76 At initialization time, a table of colors is generated. An existing color can
77 be retrieved by its index:
78 
79 ~~~ {.cpp}
80  TColor *color = gROOT->GetColor(10);
81 ~~~
82 
83 Then it can be manipulated. For example its RGB components can be modified:
84 
85 ~~~ {.cpp}
86  color->SetRGB(0.1, 0.2, 0.3);
87 ~~~
88 
89 A new color can be created the following way:
90 
91 ~~~ {.cpp}
92  Int_t ci = 1756; // color index
93  TColor *color = new TColor(ci, 0.1, 0.2, 0.3);
94 ~~~
95 
96 \since **6.07/07:**
97 TColor::GetFreeColorIndex() allows to make sure the new color is created with an
98 unused color index:
99 
100 ~~~ {.cpp}
101  Int_t ci = TColor::GetFreeColorIndex();
102  TColor *color = new TColor(ci, 0.1, 0.2, 0.3);
103 ~~~
104 
105 Two sets of colors are initialized;
106 
107  - The basic colors: colors with index from 0 to 50.
108  - The color wheel: colors with indices from 300 to 1000.
109 
110 ## <a name="C01"></a> Basic colors
111 The following image displays the 50 basic colors.
112 
113 Begin_Macro(source)
114 {
115  TCanvas *c = new TCanvas("c","Fill Area colors",0,0,500,200);
116  c->DrawColorTable();
117  return c;
118 }
119 End_Macro
120 
121 ## <a name="C02"></a> The color wheel
122 The wheel contains the recommended 216 colors to be used in web applications.
123 
124 The colors in the color wheel are created by `TColor::CreateColorWheel`.
125 
126 Using this color set for your text, background or graphics will give your
127 application a consistent appearance across different platforms and browsers.
128 
129 Colors are grouped by hue, the aspect most important in human perception.
130 Touching color chips have the same hue, but with different brightness and
131 vividness.
132 
133 Colors of slightly different hues clash. If you intend to display
134 colors of the same hue together, you should pick them from the same group.
135 
136 Each color chip is identified by a mnemonic (e.g. kYellow) and a number.
137 The keywords, kRed, kBlue, kYellow, kPink, etc are defined in the header file
138 Rtypes.h that is included in all ROOT other header files. It is better
139 to use these keywords in user code instead of hardcoded color numbers, e.g.:
140 
141 ~~~ {.cpp}
142  myObject.SetFillColor(kRed);
143  myObject.SetFillColor(kYellow-10);
144  myLine.SetLineColor(kMagenta+2);
145 ~~~
146 
147 Begin_Macro(source)
148 {
149  TColorWheel *w = new TColorWheel();
150  cw = new TCanvas("cw","cw",0,0,400,400);
151  w->SetCanvas(cw);
152  w->Draw();
153 }
154 End_Macro
155 
156 The complete list of predefined color names is the following:
157 
158 ~~~ {.cpp}
159 kWhite = 0, kBlack = 1, kGray = 920, kRed = 632, kGreen = 416,
160 kBlue = 600, kYellow = 400, kMagenta = 616, kCyan = 432, kOrange = 800,
161 kSpring = 820, kTeal = 840, kAzure = 860, kViolet = 880, kPink = 900
162 ~~~
163 
164 Note the special role of color `kWhite` (color number 0). It is the default
165 background color also. For instance in a PDF or PS files (as paper is usually white)
166 it is simply not painted. To have a white color behaving like the other color the
167 simplest is to define an other white color not attached to the color index 0:
168 
169 ~~~ {.cpp}
170  Int_t ci = TColor::GetFreeColorIndex();
171  TColor *color = new TColor(ci, 1., 1., 1.);
172 ~~~
173 
174 ## <a name="C03"></a> Bright and dark colors
175 The dark and bright color are used to give 3-D effects when drawing various
176 boxes (see TWbox, TPave, TPaveText, TPaveLabel, etc).
177 
178  - The dark colors have an index = color_index+100
179  - The bright colors have an index = color_index+150
180  - Two static functions return the bright and dark color number
181  corresponding to a color index. If the bright or dark color does not
182  exist, they are created:
183  ~~~ {.cpp}
184  Int_t dark = TColor::GetColorDark(color_index);
185  Int_t bright = TColor::GetColorBright(color_index);
186  ~~~
187 
188 ## <a name="C04"></a> Grayscale view of of canvas with colors
189 One can toggle between a grayscale preview and the regular colored mode using
190 `TCanvas::SetGrayscale()`. Note that in grayscale mode, access via RGB
191 will return grayscale values according to ITU standards (and close to b&w
192 printer gray-scales), while access via HLS returns de-saturated gray-scales. The
193 image below shows the ROOT color wheel in grayscale mode.
194 
195 Begin_Macro(source)
196 {
197  TColorWheel *w = new TColorWheel();
198  cw = new TCanvas("cw","cw",0,0,400,400);
199  cw->GetCanvas()->SetGrayscale();
200  w->SetCanvas(cw);
201  w->Draw();
202 }
203 End_Macro
204 
205 ## <a name="C05"></a> Color palettes
206 It is often very useful to represent a variable with a color map. The concept
207 of "color palette" allows to do that. One color palette is active at any time.
208 This "current palette" is set using:
209 
210 ~~~ {.cpp}
211 gStyle->SetPalette(...);
212 ~~~
213 
214 This function has two parameters: the number of colors in the palette and an
215 array of containing the indices of colors in the palette. The following small
216 example demonstrates how to define and use the color palette:
217 
218 Begin_Macro(source)
219 {
220  TCanvas *c1 = new TCanvas("c1","c1",0,0,600,400);
221  TF2 *f1 = new TF2("f1","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",1,3,1,3);
222  Int_t palette[5];
223  palette[0] = 15;
224  palette[1] = 20;
225  palette[2] = 23;
226  palette[3] = 30;
227  palette[4] = 32;
228  gStyle->SetPalette(5,palette);
229  f1->Draw("colz");
230  return c1;
231 }
232 End_Macro
233 
234  To define more a complex palette with a continuous gradient of color, one
235 should use the static function `TColor::CreateGradientColorTable()`.
236 The following example demonstrates how to proceed:
237 
238 Begin_Macro(source)
239 {
240  TCanvas *c2 = new TCanvas("c2","c2",0,0,600,400);
241  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",1,3,1,3);
242  const Int_t Number = 3;
243  Double_t Red[Number] = { 1.00, 0.00, 0.00};
244  Double_t Green[Number] = { 0.00, 1.00, 0.00};
245  Double_t Blue[Number] = { 1.00, 0.00, 1.00};
246  Double_t Length[Number] = { 0.00, 0.50, 1.00 };
247  Int_t nb=50;
248  TColor::CreateGradientColorTable(Number,Length,Red,Green,Blue,nb);
249  f2->SetContour(nb);
250  f2->SetLineWidth(1);
251  f2->SetLineColor(kBlack);
252  f2->Draw("surf1z");
253  return c2;
254 }
255 End_Macro
256 
257 The function `TColor::CreateGradientColorTable()` automatically
258 calls `gStyle->SetPalette()`, so there is not need to add one.
259 
260 After a call to `TColor::CreateGradientColorTable()` it is sometimes
261 useful to store the newly create palette for further use. In particular, it is
262 recommended to do if one wants to switch between several user define palettes.
263 To store a palette in an array it is enough to do:
264 
265 ~~~ {.cpp}
266  Int_t MyPalette[100];
267  Double_t Red[] = {0., 0.0, 1.0, 1.0, 1.0};
268  Double_t Green[] = {0., 0.0, 0.0, 1.0, 1.0};
269  Double_t Blue[] = {0., 1.0, 0.0, 0.0, 1.0};
270  Double_t Length[] = {0., .25, .50, .75, 1.0};
271  Int_t FI = TColor::CreateGradientColorTable(5, Length, Red, Green, Blue, 100);
272  for (int i=0;i<100;i++) MyPalette[i] = FI+i;
273 ~~~
274 
275 Later on to reuse the palette `MyPalette` it will be enough to do
276 
277 ~~~ {.cpp}
278  gStyle->SetPalette(100, MyPalette);
279 ~~~
280 
281 As only one palette is active, one need to use `TExec` to be able to
282 display plots using different palettes on the same pad.
283 The tutorial multipalette.C illustrates this feature.
284 
285 Begin_Macro(source)
286 ../../../tutorials/graphs/multipalette.C
287 End_Macro
288 
289 ## <a name="C06"></a> High quality predefined palettes
290 \since **6.04:**
291 62 high quality palettes are predefined with 255 colors each.
292 Despite the [disadvantages of the Rainbow color map](https://root.cern.ch/rainbow-color-map),
293 it was kept in the list of predefined color maps.
294 These palettes can be accessed "by name" with `gStyle->SetPalette(num)`.
295 `num` can be taken within the following enum:
296 
297 ~~~ {.cpp}
298 kDeepSea=51, kGreyScale=52, kDarkBodyRadiator=53,
299 kBlueYellow= 54, kRainBow=55, kInvertedDarkBodyRadiator=56,
300 kBird=57, kCubehelix=58, kGreenRedViolet=59,
301 kBlueRedYellow=60, kOcean=61, kColorPrintableOnGrey=62,
302 kAlpine=63, kAquamarine=64, kArmy=65,
303 kAtlantic=66, kAurora=67, kAvocado=68,
304 kBeach=69, kBlackBody=70, kBlueGreenYellow=71,
305 kBrownCyan=72, kCMYK=73, kCandy=74,
306 kCherry=75, kCoffee=76, kDarkRainBow=77,
307 kDarkTerrain=78, kFall=79, kFruitPunch=80,
308 kFuchsia=81, kGreyYellow=82, kGreenBrownTerrain=83,
309 kGreenPink=84, kIsland=85, kLake=86,
310 kLightTemperature=87, kLightTerrain=88, kMint=89,
311 kNeon=90, kPastel=91, kPearl=92,
312 kPigeon=93, kPlum=94, kRedBlue=95,
313 kRose=96, kRust=97, kSandyTerrain=98,
314 kSienna=99, kSolar=100, kSouthWest=101,
315 kStarryNight=102, kSunset=103, kTemperatureMap=104,
316 kThermometer=105, kValentine=106, kVisibleSpectrum=107,
317 kWaterMelon=108, kCool=109, kCopper=110,
318 kGistEarth=111, kViridis=112, kCividis=113
319 ~~~
320 
321 <table border=0>
322 <tr><td>
323 Begin_Macro
324 {
325  c = new TCanvas("c","c",0,0,300,300);
326  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
327  f2->SetContour(99); gStyle->SetPalette(kDeepSea);
328  f2->Draw("surf2Z"); f2->SetTitle("kDeepSea");
329 }
330 End_Macro
331 </td><td>
332 Begin_Macro
333 {
334  c = new TCanvas("c","c",0,0,300,300);
335  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
336  f2->SetContour(99); gStyle->SetPalette(kGreyScale);
337  f2->Draw("surf2Z"); f2->SetTitle("kGreyScale");
338 }
339 End_Macro
340 </td><td>
341 Begin_Macro
342 {
343  c = new TCanvas("c","c",0,0,300,300);
344  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
345  f2->SetContour(99); gStyle->SetPalette(kDarkBodyRadiator);
346  f2->Draw("surf2Z"); f2->SetTitle("kDarkBodyRadiator");
347 }
348 End_Macro
349 </td></tr>
350 <tr><td>
351 Begin_Macro
352 {
353  c = new TCanvas("c","c",0,0,300,300);
354  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
355  f2->SetContour(99); gStyle->SetPalette(kBlueYellow);
356  f2->Draw("surf2Z"); f2->SetTitle("kBlueYellow");
357 }
358 End_Macro
359 </td><td>
360 Begin_Macro
361 {
362  c = new TCanvas("c","c",0,0,300,300);
363  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
364  f2->SetContour(99); gStyle->SetPalette(kRainBow);
365  f2->Draw("surf2Z"); f2->SetTitle("kRainBow");
366 }
367 End_Macro
368 </td><td>
369 Begin_Macro
370 {
371  c = new TCanvas("c","c",0,0,300,300);
372  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
373  f2->SetContour(99); gStyle->SetPalette(kInvertedDarkBodyRadiator);
374  f2->Draw("surf2Z"); f2->SetTitle("kInvertedDarkBodyRadiator");
375 }
376 End_Macro
377 </td></tr>
378 <tr><td>
379 Begin_Macro
380 {
381  c = new TCanvas("c","c",0,0,300,300);
382  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
383  f2->SetContour(99); gStyle->SetPalette(kBird);
384  f2->Draw("surf2Z"); f2->SetTitle("kBird (default)");
385 }
386 End_Macro
387 </td><td>
388 Begin_Macro
389 {
390  c = new TCanvas("c","c",0,0,300,300);
391  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
392  f2->SetContour(99); gStyle->SetPalette(kCubehelix);
393  f2->Draw("surf2Z"); f2->SetTitle("kCubehelix");
394 }
395 End_Macro
396 </td><td>
397 Begin_Macro
398 {
399  c = new TCanvas("c","c",0,0,300,300);
400  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
401  f2->SetContour(99); gStyle->SetPalette(kGreenRedViolet);
402  f2->Draw("surf2Z"); f2->SetTitle("kGreenRedViolet");
403 }
404 End_Macro
405 </td></tr>
406 <tr><td>
407 Begin_Macro
408 {
409  c = new TCanvas("c","c",0,0,300,300);
410  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
411  f2->SetContour(99); gStyle->SetPalette(kBlueRedYellow);
412  f2->Draw("surf2Z"); f2->SetTitle("kBlueRedYellow");
413 }
414 End_Macro
415 </td><td>
416 Begin_Macro
417 {
418  c = new TCanvas("c","c",0,0,300,300);
419  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
420  f2->SetContour(99); gStyle->SetPalette(kOcean);
421  f2->Draw("surf2Z"); f2->SetTitle("kOcean");
422 }
423 End_Macro
424 </td><td>
425 Begin_Macro
426 {
427  c = new TCanvas("c","c",0,0,300,300);
428  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
429  f2->SetContour(99); gStyle->SetPalette(kColorPrintableOnGrey);
430  f2->Draw("surf2Z"); f2->SetTitle("kColorPrintableOnGrey");
431 }
432 End_Macro
433 </td></tr>
434 <tr><td>
435 Begin_Macro
436 {
437  c = new TCanvas("c","c",0,0,300,300);
438  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
439  f2->SetContour(99); gStyle->SetPalette(kAlpine);
440  f2->Draw("surf2Z"); f2->SetTitle("kAlpine");
441 }
442 End_Macro
443 </td><td>
444 Begin_Macro
445 {
446  c = new TCanvas("c","c",0,0,300,300);
447  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
448  f2->SetContour(99); gStyle->SetPalette(kAquamarine);
449  f2->Draw("surf2Z"); f2->SetTitle("kAquamarine");
450 }
451 End_Macro
452 </td><td>
453 Begin_Macro
454 {
455  c = new TCanvas("c","c",0,0,300,300);
456  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
457  f2->SetContour(99); gStyle->SetPalette(kArmy);
458  f2->Draw("surf2Z"); f2->SetTitle("kArmy");
459 }
460 End_Macro
461 </td></tr>
462 <tr><td>
463 Begin_Macro
464 {
465  c = new TCanvas("c","c",0,0,300,300);
466  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
467  f2->SetContour(99); gStyle->SetPalette(kAtlantic);
468  f2->Draw("surf2Z"); f2->SetTitle("kAtlantic");
469 }
470 End_Macro
471 </td><td>
472 Begin_Macro
473 {
474  c = new TCanvas("c","c",0,0,300,300);
475  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
476  f2->SetContour(99); gStyle->SetPalette(kAurora);
477  f2->Draw("surf2Z"); f2->SetTitle("kAurora");
478 }
479 End_Macro
480 </td><td>
481 Begin_Macro
482 {
483  c = new TCanvas("c","c",0,0,300,300);
484  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
485  f2->SetContour(99); gStyle->SetPalette(kAvocado);
486  f2->Draw("surf2Z"); f2->SetTitle("kAvocado");
487 }
488 End_Macro
489 </td></tr>
490 <tr><td>
491 Begin_Macro
492 {
493  c = new TCanvas("c","c",0,0,300,300);
494  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
495  f2->SetContour(99); gStyle->SetPalette(kBeach);
496  f2->Draw("surf2Z"); f2->SetTitle("kBeach");
497 }
498 End_Macro
499 </td><td>
500 Begin_Macro
501 {
502  c = new TCanvas("c","c",0,0,300,300);
503  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
504  f2->SetContour(99); gStyle->SetPalette(kBlackBody);
505  f2->Draw("surf2Z"); f2->SetTitle("kBlackBody");
506 }
507 End_Macro
508 </td><td>
509 Begin_Macro
510 {
511  c = new TCanvas("c","c",0,0,300,300);
512  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
513  f2->SetContour(99); gStyle->SetPalette(kBlueGreenYellow);
514  f2->Draw("surf2Z"); f2->SetTitle("kBlueGreenYellow");
515 }
516 End_Macro
517 </td></tr>
518 <tr><td>
519 Begin_Macro
520 {
521  c = new TCanvas("c","c",0,0,300,300);
522  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
523  f2->SetContour(99); gStyle->SetPalette(kBrownCyan);
524  f2->Draw("surf2Z"); f2->SetTitle("kBrownCyan");
525 }
526 End_Macro
527 </td><td>
528 Begin_Macro
529 {
530  c = new TCanvas("c","c",0,0,300,300);
531  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
532  f2->SetContour(99); gStyle->SetPalette(kCMYK);
533  f2->Draw("surf2Z"); f2->SetTitle("kCMYK");
534 }
535 End_Macro
536 </td><td>
537 Begin_Macro
538 {
539  c = new TCanvas("c","c",0,0,300,300);
540  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
541  f2->SetContour(99); gStyle->SetPalette(kCandy);
542  f2->Draw("surf2Z"); f2->SetTitle("kCandy");
543 }
544 End_Macro
545 </td></tr>
546 <tr><td>
547 Begin_Macro
548 {
549  c = new TCanvas("c","c",0,0,300,300);
550  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
551  f2->SetContour(99); gStyle->SetPalette(kCherry);
552  f2->Draw("surf2Z"); f2->SetTitle("kCherry");
553 }
554 End_Macro
555 </td><td>
556 Begin_Macro
557 {
558  c = new TCanvas("c","c",0,0,300,300);
559  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
560  f2->SetContour(99); gStyle->SetPalette(kCoffee);
561  f2->Draw("surf2Z"); f2->SetTitle("kCoffee");
562 }
563 End_Macro
564 </td><td>
565 Begin_Macro
566 {
567  c = new TCanvas("c","c",0,0,300,300);
568  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
569  f2->SetContour(99); gStyle->SetPalette(kDarkRainBow);
570  f2->Draw("surf2Z"); f2->SetTitle("kDarkRainBow");
571 }
572 End_Macro
573 </td></tr>
574 <tr><td>
575 Begin_Macro
576 {
577  c = new TCanvas("c","c",0,0,300,300);
578  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
579  f2->SetContour(99); gStyle->SetPalette(kDarkTerrain);
580  f2->Draw("surf2Z"); f2->SetTitle("kDarkTerrain");
581 }
582 End_Macro
583 </td><td>
584 Begin_Macro
585 {
586  c = new TCanvas("c","c",0,0,300,300);
587  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
588  f2->SetContour(99); gStyle->SetPalette(kFall);
589  f2->Draw("surf2Z"); f2->SetTitle("kFall");
590 }
591 End_Macro
592 </td><td>
593 Begin_Macro
594 {
595  c = new TCanvas("c","c",0,0,300,300);
596  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
597  f2->SetContour(99); gStyle->SetPalette(kFruitPunch);
598  f2->Draw("surf2Z"); f2->SetTitle("kFruitPunch");
599 }
600 End_Macro
601 </td></tr>
602 <tr><td>
603 Begin_Macro
604 {
605  c = new TCanvas("c","c",0,0,300,300);
606  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
607  f2->SetContour(99); gStyle->SetPalette(kFuchsia);
608  f2->Draw("surf2Z"); f2->SetTitle("kFuchsia");
609 }
610 End_Macro
611 </td><td>
612 Begin_Macro
613 {
614  c = new TCanvas("c","c",0,0,300,300);
615  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
616  f2->SetContour(99); gStyle->SetPalette(kGreyYellow);
617  f2->Draw("surf2Z"); f2->SetTitle("kGreyYellow");
618 }
619 End_Macro
620 </td><td>
621 Begin_Macro
622 {
623  c = new TCanvas("c","c",0,0,300,300);
624  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
625  f2->SetContour(99); gStyle->SetPalette(kGreenBrownTerrain);
626  f2->Draw("surf2Z"); f2->SetTitle("kGreenBrownTerrain");
627 }
628 End_Macro
629 </td></tr>
630 <tr><td>
631 Begin_Macro
632 {
633  c = new TCanvas("c","c",0,0,300,300);
634  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
635  f2->SetContour(99); gStyle->SetPalette(kGreenPink);
636  f2->Draw("surf2Z"); f2->SetTitle("kGreenPink");
637 }
638 End_Macro
639 </td><td>
640 Begin_Macro
641 {
642  c = new TCanvas("c","c",0,0,300,300);
643  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
644  f2->SetContour(99); gStyle->SetPalette(kIsland);
645  f2->Draw("surf2Z"); f2->SetTitle("kIsland");
646 }
647 End_Macro
648 </td><td>
649 Begin_Macro
650 {
651  c = new TCanvas("c","c",0,0,300,300);
652  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
653  f2->SetContour(99); gStyle->SetPalette(kLake);
654  f2->Draw("surf2Z"); f2->SetTitle("kLake");
655 }
656 End_Macro
657 </td></tr>
658 <tr><td>
659 Begin_Macro
660 {
661  c = new TCanvas("c","c",0,0,300,300);
662  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
663  f2->SetContour(99); gStyle->SetPalette(kLightTemperature);
664  f2->Draw("surf2Z"); f2->SetTitle("kLightTemperature");
665 }
666 End_Macro
667 </td><td>
668 Begin_Macro
669 {
670  c = new TCanvas("c","c",0,0,300,300);
671  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
672  f2->SetContour(99); gStyle->SetPalette(kLightTerrain);
673  f2->Draw("surf2Z"); f2->SetTitle("kLightTerrain");
674 }
675 End_Macro
676 </td><td>
677 Begin_Macro
678 {
679  c = new TCanvas("c","c",0,0,300,300);
680  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
681  f2->SetContour(99); gStyle->SetPalette(kMint);
682  f2->Draw("surf2Z"); f2->SetTitle("kMint");
683 }
684 End_Macro
685 </td></tr>
686 <tr><td>
687 Begin_Macro
688 {
689  c = new TCanvas("c","c",0,0,300,300);
690  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
691  f2->SetContour(99); gStyle->SetPalette(kNeon);
692  f2->Draw("surf2Z"); f2->SetTitle("kNeon");
693 }
694 End_Macro
695 </td><td>
696 Begin_Macro
697 {
698  c = new TCanvas("c","c",0,0,300,300);
699  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
700  f2->SetContour(99); gStyle->SetPalette(kPastel);
701  f2->Draw("surf2Z"); f2->SetTitle("kPastel");
702 }
703 End_Macro
704 </td><td>
705 Begin_Macro
706 {
707  c = new TCanvas("c","c",0,0,300,300);
708  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
709  f2->SetContour(99); gStyle->SetPalette(kPearl);
710  f2->Draw("surf2Z"); f2->SetTitle("kPearl");
711 }
712 End_Macro
713 </td></tr>
714 <tr><td>
715 Begin_Macro
716 {
717  c = new TCanvas("c","c",0,0,300,300);
718  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
719  f2->SetContour(99); gStyle->SetPalette(kPigeon);
720  f2->Draw("surf2Z"); f2->SetTitle("kPigeon");
721 }
722 End_Macro
723 </td><td>
724 Begin_Macro
725 {
726  c = new TCanvas("c","c",0,0,300,300);
727  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
728  f2->SetContour(99); gStyle->SetPalette(kPlum);
729  f2->Draw("surf2Z"); f2->SetTitle("kPlum");
730 }
731 End_Macro
732 </td><td>
733 Begin_Macro
734 {
735  c = new TCanvas("c","c",0,0,300,300);
736  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
737  f2->SetContour(99); gStyle->SetPalette(kRedBlue);
738  f2->Draw("surf2Z"); f2->SetTitle("kRedBlue");
739 }
740 End_Macro
741 </td></tr>
742 <tr><td>
743 Begin_Macro
744 {
745  c = new TCanvas("c","c",0,0,300,300);
746  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
747  f2->SetContour(99); gStyle->SetPalette(kRose);
748  f2->Draw("surf2Z"); f2->SetTitle("kRose");
749 }
750 End_Macro
751 </td><td>
752 Begin_Macro
753 {
754  c = new TCanvas("c","c",0,0,300,300);
755  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
756  f2->SetContour(99); gStyle->SetPalette(kRust);
757  f2->Draw("surf2Z"); f2->SetTitle("kRust");
758 }
759 End_Macro
760 </td><td>
761 Begin_Macro
762 {
763  c = new TCanvas("c","c",0,0,300,300);
764  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
765  f2->SetContour(99); gStyle->SetPalette(kSandyTerrain);
766  f2->Draw("surf2Z"); f2->SetTitle("kSandyTerrain");
767 }
768 End_Macro
769 </td></tr>
770 <tr><td>
771 Begin_Macro
772 {
773  c = new TCanvas("c","c",0,0,300,300);
774  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
775  f2->SetContour(99); gStyle->SetPalette(kSienna);
776  f2->Draw("surf2Z"); f2->SetTitle("kSienna");
777 }
778 End_Macro
779 </td><td>
780 Begin_Macro
781 {
782  c = new TCanvas("c","c",0,0,300,300);
783  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
784  f2->SetContour(99); gStyle->SetPalette(kSolar);
785  f2->Draw("surf2Z"); f2->SetTitle("kSolar");
786 }
787 End_Macro
788 </td><td>
789 Begin_Macro
790 {
791  c = new TCanvas("c","c",0,0,300,300);
792  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
793  f2->SetContour(99); gStyle->SetPalette(kSouthWest);
794  f2->Draw("surf2Z"); f2->SetTitle("kSouthWest");
795 }
796 End_Macro
797 </td></tr>
798 <tr><td>
799 Begin_Macro
800 {
801  c = new TCanvas("c","c",0,0,300,300);
802  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
803  f2->SetContour(99); gStyle->SetPalette(kStarryNight);
804  f2->Draw("surf2Z"); f2->SetTitle("kStarryNight");
805 }
806 End_Macro
807 </td><td>
808 Begin_Macro
809 {
810  c = new TCanvas("c","c",0,0,300,300);
811  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
812  f2->SetContour(99); gStyle->SetPalette(kSunset);
813  f2->Draw("surf2Z"); f2->SetTitle("kSunset");
814 }
815 End_Macro
816 </td><td>
817 Begin_Macro
818 {
819  c = new TCanvas("c","c",0,0,300,300);
820  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
821  f2->SetContour(99); gStyle->SetPalette(kTemperatureMap);
822  f2->Draw("surf2Z"); f2->SetTitle("kTemperatureMap");
823 }
824 End_Macro
825 </td></tr>
826 <tr><td>
827 Begin_Macro
828 {
829  c = new TCanvas("c","c",0,0,300,300);
830  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
831  f2->SetContour(99); gStyle->SetPalette(kThermometer);
832  f2->Draw("surf2Z"); f2->SetTitle("kThermometer");
833 }
834 End_Macro
835 </td><td>
836 Begin_Macro
837 {
838  c = new TCanvas("c","c",0,0,300,300);
839  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
840  f2->SetContour(99); gStyle->SetPalette(kValentine);
841  f2->Draw("surf2Z"); f2->SetTitle("kValentine");
842 }
843 End_Macro
844 </td><td>
845 Begin_Macro
846 {
847  c = new TCanvas("c","c",0,0,300,300);
848  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
849  f2->SetContour(99); gStyle->SetPalette(kVisibleSpectrum);
850  f2->Draw("surf2Z"); f2->SetTitle("kVisibleSpectrum");
851 }
852 End_Macro
853 </td></tr>
854 <tr><td>
855 Begin_Macro
856 {
857  c = new TCanvas("c","c",0,0,300,300);
858  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
859  f2->SetContour(99); gStyle->SetPalette(kWaterMelon);
860  f2->Draw("surf2Z"); f2->SetTitle("kWaterMelon");
861 }
862 End_Macro
863 </td><td>
864 Begin_Macro
865 {
866  c = new TCanvas("c","c",0,0,300,300);
867  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
868  f2->SetContour(99); gStyle->SetPalette(kCool);
869  f2->Draw("surf2Z"); f2->SetTitle("kCool");
870 }
871 End_Macro
872 </td><td>
873 Begin_Macro
874 {
875  c = new TCanvas("c","c",0,0,300,300);
876  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
877  f2->SetContour(99); gStyle->SetPalette(kCopper);
878  f2->Draw("surf2Z"); f2->SetTitle("kCopper");
879 }
880 End_Macro
881 </td></tr>
882 <tr><td>
883 Begin_Macro
884 {
885  c = new TCanvas("c","c",0,0,300,300);
886  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
887  f2->SetContour(99); gStyle->SetPalette(kGistEarth);
888  f2->Draw("surf2Z"); f2->SetTitle("kGistEarth");
889 }
890 End_Macro
891 </td><td>
892 Begin_Macro
893 {
894  c = new TCanvas("c","c",0,0,300,300);
895  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
896  f2->SetContour(99); gStyle->SetPalette(kViridis);
897  f2->Draw("surf2Z"); f2->SetTitle("kViridis");
898 }
899 End_Macro
900 </td><td>
901 Begin_Macro
902 {
903  c = new TCanvas("c","c",0,0,300,300);
904  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
905  f2->SetContour(99); gStyle->SetPalette(kCividis);
906  f2->Draw("surf2Z"); f2->SetTitle("kCividis");
907 }
908 End_Macro
909 </td></tr>
910 </table>
911 
912 ## <a name="C061"></a> Palette inversion
913 Once a palette is defined, it is possible to invert the color order thanks to the
914 method TColor::InvertPalette. The top of the palette becomes the bottom and vice versa.
915 
916 Begin_Macro(source)
917 {
918  auto c = new TCanvas("c","c",0,0,600,400);
919  TF2 *f2 = new TF2("f2","0.1+(1-(x-2)*(x-2))*(1-(y-2)*(y-2))",0.999,3.002,0.999,3.002);
920  f2->SetContour(99); gStyle->SetPalette(kCherry);
921  TColor::InvertPalette();
922  f2->Draw("surf2Z"); f2->SetTitle("kCherry inverted");
923 }
924 End_Macro
925 
926 ## <a name="C07"></a> Color transparency
927 To make a graphics object transparent it is enough to set its color to a
928 transparent one. The color transparency is defined via its alpha component. The
929 alpha value varies from `0.` (fully transparent) to `1.` (fully
930 opaque). To set the alpha value of an existing color it is enough to do:
931 
932 ~~~ {.cpp}
933  TColor *col26 = gROOT->GetColor(26);
934  col26->SetAlpha(0.01);
935 ~~~
936 
937 A new color can be created transparent the following way:
938 
939 ~~~ {.cpp}
940  Int_t ci = 1756;
941  TColor *color = new TColor(ci, 0.1, 0.2, 0.3, "", 0.5); // alpha = 0.5
942 ~~~
943 
944 An example of transparency usage with parallel coordinates can be found
945 in parallelcoordtrans.C.
946 
947 To ease the creation of a transparent color the static method
948 `GetColorTransparent(Int_t color, Float_t a)` is provided.
949 In the following example the `trans_red` color index point to
950 a red color 30% transparent. The alpha value of the color index
951 `kRed` is not modified.
952 
953 ~~~ {.cpp}
954  Int_t trans_red = GetColorTransparent(kRed, 0.3);
955 ~~~
956 
957 This function is also used in the methods
958 `SetFillColorAlpha()`, `SetLineColorAlpha()`,
959 `SetMarkerColorAlpha()` and `SetTextColorAlpha()`.
960 In the following example the fill color of the histogram `histo`
961 is set to blue with a transparency of 35%. The color `kBlue`
962 itself remains fully opaque.
963 
964 ~~~ {.cpp}
965  histo->SetFillColorAlpha(kBlue, 0.35);
966 ~~~
967 
968 The transparency is available on all platforms when the flag `OpenGL.CanvasPreferGL` is set to `1`
969 in `$ROOTSYS/etc/system.rootrc`, or on Mac with the Cocoa backend. On the file output
970 it is visible with PDF, PNG, Gif, JPEG, SVG, TeX ... but not PostScript.
971 The following macro gives an example of transparency usage:
972 
973 Begin_Macro(source)
974 ../../../tutorials/graphics/transparency.C
975 End_Macro
976 
977 */
978 
979 ////////////////////////////////////////////////////////////////////////////////
980 /// Default constructor.
981 
983 {
984  fNumber = -1;
985  fRed = fGreen = fBlue = fHue = fLight = fSaturation = -1;
986  fAlpha = 1;
987 }
988 
989 ////////////////////////////////////////////////////////////////////////////////
990 /// Normal color constructor. Initialize a color structure.
991 /// Compute the RGB and HLS color components.
992 
994  Float_t a)
995  : TNamed(name,"")
996 {
998  // do not enter if color number already exist
999  TColor *col = gROOT->GetColor(color);
1000  if (col) {
1001  Warning("TColor", "color %d already defined", color);
1002  fNumber = col->GetNumber();
1003  fRed = col->GetRed();
1004  fGreen = col->GetGreen();
1005  fBlue = col->GetBlue();
1006  fHue = col->GetHue();
1007  fLight = col->GetLight();
1008  fAlpha = col->GetAlpha();
1009  fSaturation = col->GetSaturation();
1010  return;
1011  }
1012 
1013  fNumber = color;
1014 
1016 
1017  char aname[32];
1018  if (!name || !*name) {
1019  snprintf(aname,32, "Color%d", color);
1020  SetName(aname);
1021  }
1022 
1023  // enter in the list of colors
1024  TObjArray *lcolors = (TObjArray*)gROOT->GetListOfColors();
1025  lcolors->AddAtAndExpand(this, color);
1026 
1027  // fill color structure
1028  SetRGB(r, g, b);
1029  fAlpha = a;
1030  gDefinedColors++;
1031 }
1032 
1033 ////////////////////////////////////////////////////////////////////////////////
1034 /// Fast TColor constructor. It creates a color with an index just above the
1035 /// current highest one. It does not name the color.
1036 /// This is useful to create palettes.
1037 
1039 {
1042  fRed = r;
1043  fGreen = g;
1044  fBlue = b;
1045  fAlpha = a;
1047 
1048  // enter in the list of colors
1049  TObjArray *lcolors = (TObjArray*)gROOT->GetListOfColors();
1050  lcolors->AddAtAndExpand(this, fNumber);
1051  gDefinedColors++;
1052 }
1053 
1054 ////////////////////////////////////////////////////////////////////////////////
1055 /// Color destructor.
1056 
1058 {
1059  gROOT->GetListOfColors()->Remove(this);
1060  if (gROOT->GetListOfColors()->GetEntries() == 0) {fgPalette.Set(0); fgPalette=0;}
1061 }
1062 
1063 ////////////////////////////////////////////////////////////////////////////////
1064 /// Color copy constructor.
1065 
1066 TColor::TColor(const TColor &color) : TNamed(color)
1067 {
1068  ((TColor&)color).Copy(*this);
1069 }
1070 
1072 {
1073  ((TColor &)color).Copy(*this);
1074  return *this;
1075 }
1076 
1077 ////////////////////////////////////////////////////////////////////////////////
1078 /// Initialize colors used by the TCanvas based graphics (via TColor objects).
1079 /// This method should be called before the ApplicationImp is created (which
1080 /// initializes the GUI colors).
1081 
1083 {
1084  static Bool_t initDone = kFALSE;
1085 
1086  if (initDone) return;
1087  initDone = kTRUE;
1088 
1089  if (gROOT->GetListOfColors()->First() == 0) {
1090 
1091  new TColor(kWhite,1,1,1,"background");
1092  new TColor(kBlack,0,0,0,"black");
1093  new TColor(2,1,0,0,"red");
1094  new TColor(3,0,1,0,"green");
1095  new TColor(4,0,0,1,"blue");
1096  new TColor(5,1,1,0,"yellow");
1097  new TColor(6,1,0,1,"magenta");
1098  new TColor(7,0,1,1,"cyan");
1099  new TColor(10,0.999,0.999,0.999,"white");
1100  new TColor(11,0.754,0.715,0.676,"editcol");
1101 
1102  // The color white above is defined as being nearly white.
1103  // Sets the associated dark color also to white.
1105  TColor *c110 = gROOT->GetColor(110);
1106  if (c110) c110->SetRGB(0.999,0.999,.999);
1107 
1108  // Initialize Custom colors
1109  new TColor(20,0.8,0.78,0.67);
1110  new TColor(31,0.54,0.66,0.63);
1111  new TColor(41,0.83,0.81,0.53);
1112  new TColor(30,0.52,0.76,0.64);
1113  new TColor(32,0.51,0.62,0.55);
1114  new TColor(24,0.70,0.65,0.59);
1115  new TColor(21,0.8,0.78,0.67);
1116  new TColor(47,0.67,0.56,0.58);
1117  new TColor(35,0.46,0.54,0.57);
1118  new TColor(33,0.68,0.74,0.78);
1119  new TColor(39,0.5,0.5,0.61);
1120  new TColor(37,0.43,0.48,0.52);
1121  new TColor(38,0.49,0.6,0.82);
1122  new TColor(36,0.41,0.51,0.59);
1123  new TColor(49,0.58,0.41,0.44);
1124  new TColor(43,0.74,0.62,0.51);
1125  new TColor(22,0.76,0.75,0.66);
1126  new TColor(45,0.75,0.51,0.47);
1127  new TColor(44,0.78,0.6,0.49);
1128  new TColor(26,0.68,0.6,0.55);
1129  new TColor(28,0.53,0.4,0.34);
1130  new TColor(25,0.72,0.64,0.61);
1131  new TColor(27,0.61,0.56,0.51);
1132  new TColor(23,0.73,0.71,0.64);
1133  new TColor(42,0.87,0.73,0.53);
1134  new TColor(46,0.81,0.37,0.38);
1135  new TColor(48,0.65,0.47,0.48);
1136  new TColor(34,0.48,0.56,0.6);
1137  new TColor(40,0.67,0.65,0.75);
1138  new TColor(29,0.69,0.81,0.78);
1139 
1140  // Initialize some additional greyish non saturated colors
1141  new TColor(8, 0.35,0.83,0.33);
1142  new TColor(9, 0.35,0.33,0.85);
1143  new TColor(12,.3,.3,.3,"grey12");
1144  new TColor(13,.4,.4,.4,"grey13");
1145  new TColor(14,.5,.5,.5,"grey14");
1146  new TColor(15,.6,.6,.6,"grey15");
1147  new TColor(16,.7,.7,.7,"grey16");
1148  new TColor(17,.8,.8,.8,"grey17");
1149  new TColor(18,.9,.9,.9,"grey18");
1150  new TColor(19,.95,.95,.95,"grey19");
1151  new TColor(50, 0.83,0.35,0.33);
1152 
1153  // Initialize the Pretty Palette Spectrum Violet->Red
1154  // The color model used here is based on the HLS model which
1155  // is much more suitable for creating palettes than RGB.
1156  // Fixing the saturation and lightness we can scan through the
1157  // spectrum of visible light by using "hue" alone.
1158  // In Root hue takes values from 0 to 360.
1159  Int_t i;
1160  Float_t saturation = 1;
1161  Float_t lightness = 0.5;
1162  Float_t maxHue = 280;
1163  Float_t minHue = 0;
1164  Int_t maxPretty = 50;
1165  Float_t hue;
1166  Float_t r=0., g=0., b=0., h, l, s;
1167 
1168  for (i=0 ; i<maxPretty-1 ; i++) {
1169  hue = maxHue-(i+1)*((maxHue-minHue)/maxPretty);
1170  TColor::HLStoRGB(hue, lightness, saturation, r, g, b);
1171  new TColor(i+51, r, g, b);
1172  }
1173 
1174  // Initialize special colors for x3d
1175  TColor *s0;
1176  for (i = 1; i < 8; i++) {
1177  s0 = gROOT->GetColor(i);
1178  if (s0) s0->GetRGB(r,g,b);
1179  if (i == 1) { r = 0.6; g = 0.6; b = 0.6; }
1180  if (r == 1) r = 0.9; else if (r == 0) r = 0.1;
1181  if (g == 1) g = 0.9; else if (g == 0) g = 0.1;
1182  if (b == 1) b = 0.9; else if (b == 0) b = 0.1;
1183  TColor::RGBtoHLS(r,g,b,h,l,s);
1184  TColor::HLStoRGB(h,0.6*l,s,r,g,b);
1185  new TColor(200+4*i-3,r,g,b);
1186  TColor::HLStoRGB(h,0.8*l,s,r,g,b);
1187  new TColor(200+4*i-2,r,g,b);
1188  TColor::HLStoRGB(h,1.2*l,s,r,g,b);
1189  new TColor(200+4*i-1,r,g,b);
1190  TColor::HLStoRGB(h,1.4*l,s,r,g,b);
1191  new TColor(200+4*i ,r,g,b);
1192  }
1193 
1194  // Create the ROOT Color Wheel
1196  }
1197  // If fgPalette.fN !=0 SetPalette has been called already
1198  // (from rootlogon.C for instance)
1199 
1200  if (!fgPalette.fN) SetPalette(1,0);
1201 }
1202 
1203 ////////////////////////////////////////////////////////////////////////////////
1204 /// Return color as hexadecimal string. This string can be directly passed
1205 /// to, for example, TGClient::GetColorByName(). String will be reused so
1206 /// copy immediately if needed.
1207 
1208 const char *TColor::AsHexString() const
1209 {
1210  static TString tempbuf;
1211 
1212  Int_t r, g, b, a;
1213  r = Int_t(GetRed() * 255);
1214  g = Int_t(GetGreen() * 255);
1215  b = Int_t(GetBlue() * 255);
1216  a = Int_t(fAlpha * 255);
1217 
1218  if (a != 255) {
1219  tempbuf.Form("#%02x%02x%02x%02x", a, r, g, b);
1220  } else {
1221  tempbuf.Form("#%02x%02x%02x", r, g, b);
1222  }
1223  return tempbuf;
1224 }
1225 
1226 ////////////////////////////////////////////////////////////////////////////////
1227 /// Copy this color to obj.
1228 
1229 void TColor::Copy(TObject &obj) const
1230 {
1231  TNamed::Copy((TNamed&)obj);
1232  ((TColor&)obj).fRed = fRed;
1233  ((TColor&)obj).fGreen = fGreen;
1234  ((TColor&)obj).fBlue = fBlue;
1235  ((TColor&)obj).fHue = fHue;
1236  ((TColor&)obj).fLight = fLight;
1237  ((TColor&)obj).fAlpha = fAlpha;
1238  ((TColor&)obj).fSaturation = fSaturation;
1239  ((TColor&)obj).fNumber = fNumber;
1240 }
1241 
1242 ////////////////////////////////////////////////////////////////////////////////
1243 /// Create the Gray scale colors in the Color Wheel
1244 
1246 {
1247  if (gROOT->GetColor(kGray)) return;
1248  TColor *gray = new TColor(kGray,204./255.,204./255.,204./255.);
1249  TColor *gray1 = new TColor(kGray+1,153./255.,153./255.,153./255.);
1250  TColor *gray2 = new TColor(kGray+2,102./255.,102./255.,102./255.);
1251  TColor *gray3 = new TColor(kGray+3, 51./255., 51./255., 51./255.);
1252  gray ->SetName("kGray");
1253  gray1->SetName("kGray+1");
1254  gray2->SetName("kGray+2");
1255  gray3->SetName("kGray+3");
1256 }
1257 
1258 ////////////////////////////////////////////////////////////////////////////////
1259 /// Create the "circle" colors in the color wheel.
1260 
1261 void TColor::CreateColorsCircle(Int_t offset, const char *name, UChar_t *rgb)
1262 {
1263  TString colorname;
1264  for (Int_t n=0;n<15;n++) {
1265  Int_t colorn = offset+n-10;
1266  TColor *color = gROOT->GetColor(colorn);
1267  if (!color) {
1268  color = new TColor(colorn,rgb[3*n]/255.,rgb[3*n+1]/255.,rgb[3*n+2]/255.);
1269  color->SetTitle(color->AsHexString());
1270  if (n>10) colorname.Form("%s+%d",name,n-10);
1271  else if (n<10) colorname.Form("%s-%d",name,10-n);
1272  else colorname.Form("%s",name);
1273  color->SetName(colorname);
1274  }
1275  }
1276 }
1277 
1278 ////////////////////////////////////////////////////////////////////////////////
1279 /// Create the "rectangular" colors in the color wheel.
1280 
1281 void TColor::CreateColorsRectangle(Int_t offset, const char *name, UChar_t *rgb)
1282 {
1283  TString colorname;
1284  for (Int_t n=0;n<20;n++) {
1285  Int_t colorn = offset+n-9;
1286  TColor *color = gROOT->GetColor(colorn);
1287  if (!color) {
1288  color = new TColor(colorn,rgb[3*n]/255.,rgb[3*n+1]/255.,rgb[3*n+2]/255.);
1289  color->SetTitle(color->AsHexString());
1290  if (n>9) colorname.Form("%s+%d",name,n-9);
1291  else if (n<9) colorname.Form("%s-%d",name,9-n);
1292  else colorname.Form("%s",name);
1293  color->SetName(colorname);
1294  }
1295  }
1296 }
1297 
1298 ////////////////////////////////////////////////////////////////////////////////
1299 /// Static function steering the creation of all colors in the color wheel.
1300 
1302 {
1303  UChar_t magenta[46]= {255,204,255
1304  ,255,153,255, 204,153,204
1305  ,255,102,255, 204,102,204, 153,102,153
1306  ,255, 51,255, 204, 51,204, 153, 51,153, 102, 51,102
1307  ,255, 0,255, 204, 0,204, 153, 0,153, 102, 0,102, 51, 0, 51};
1308 
1309  UChar_t red[46] = {255,204,204
1310  ,255,153,153, 204,153,153
1311  ,255,102,102, 204,102,102, 153,102,102
1312  ,255, 51, 51, 204, 51, 51, 153, 51, 51, 102, 51, 51
1313  ,255, 0, 0, 204, 0, 0, 153, 0, 0, 102, 0, 0, 51, 0, 0};
1314 
1315  UChar_t yellow[46] = {255,255,204
1316  ,255,255,153, 204,204,153
1317  ,255,255,102, 204,204,102, 153,153,102
1318  ,255,255, 51, 204,204, 51, 153,153, 51, 102,102, 51
1319  ,255,255, 0, 204,204, 0, 153,153, 0, 102,102, 0, 51, 51, 0};
1320 
1321  UChar_t green[46] = {204,255,204
1322  ,153,255,153, 153,204,153
1323  ,102,255,102, 102,204,102, 102,153,102
1324  , 51,255, 51, 51,204, 51, 51,153, 51, 51,102, 51
1325  , 0,255, 0, 0,204, 0, 0,153, 0, 0,102, 0, 0, 51, 0};
1326 
1327  UChar_t cyan[46] = {204,255,255
1328  ,153,255,255, 153,204,204
1329  ,102,255,255, 102,204,204, 102,153,153
1330  , 51,255,255, 51,204,204, 51,153,153, 51,102,102
1331  , 0,255,255, 0,204,204, 0,153,153, 0,102,102, 0, 51, 51};
1332 
1333  UChar_t blue[46] = {204,204,255
1334  ,153,153,255, 153,153,204
1335  ,102,102,255, 102,102,204, 102,102,153
1336  , 51, 51,255, 51, 51,204, 51, 51,153, 51, 51,102
1337  , 0, 0,255, 0, 0,204, 0, 0,153, 0, 0,102, 0, 0, 51};
1338 
1339  UChar_t pink[60] = {255, 51,153, 204, 0,102, 102, 0, 51, 153, 0, 51, 204, 51,102
1340  ,255,102,153, 255, 0,102, 255, 51,102, 204, 0, 51, 255, 0, 51
1341  ,255,153,204, 204,102,153, 153, 51,102, 153, 0,102, 204, 51,153
1342  ,255,102,204, 255, 0,153, 204, 0,153, 255, 51,204, 255, 0,153};
1343 
1344  UChar_t orange[60]={255,204,153, 204,153,102, 153,102, 51, 153,102, 0, 204,153, 51
1345  ,255,204,102, 255,153, 0, 255,204, 51, 204,153, 0, 255,204, 0
1346  ,255,153, 51, 204,102, 0, 102, 51, 0, 153, 51, 0, 204,102, 51
1347  ,255,153,102, 255,102, 0, 255,102, 51, 204, 51, 0, 255, 51, 0};
1348 
1349  UChar_t spring[60]={153,255, 51, 102,204, 0, 51,102, 0, 51,153, 0, 102,204, 51
1350  ,153,255,102, 102,255, 0, 102,255, 51, 51,204, 0, 51,255, 0
1351  ,204,255,153, 153,204,102, 102,153, 51, 102,153, 0, 153,204, 51
1352  ,204,255,102, 153,255, 0, 204,255, 51, 153,204, 0, 204,255, 0};
1353 
1354  UChar_t teal[60] = {153,255,204, 102,204,153, 51,153,102, 0,153,102, 51,204,153
1355  ,102,255,204, 0,255,102, 51,255,204, 0,204,153, 0,255,204
1356  , 51,255,153, 0,204,102, 0,102, 51, 0,153, 51, 51,204,102
1357  ,102,255,153, 0,255,153, 51,255,102, 0,204, 51, 0,255, 51};
1358 
1359  UChar_t azure[60] ={153,204,255, 102,153,204, 51,102,153, 0, 51,153, 51,102,204
1360  ,102,153,255, 0,102,255, 51,102,255, 0, 51,204, 0, 51,255
1361  , 51,153,255, 0,102,204, 0, 51,102, 0,102,153, 51,153,204
1362  ,102,204,255, 0,153,255, 51,204,255, 0,153,204, 0,204,255};
1363 
1364  UChar_t violet[60]={204,153,255, 153,102,204, 102, 51,153, 102, 0,153, 153, 51,204
1365  ,204,102,255, 153, 0,255, 204, 51,255, 153, 0,204, 204, 0,255
1366  ,153, 51,255, 102, 0,204, 51, 0,102, 51, 0,153, 102, 51,204
1367  ,153,102,255, 102, 0,255, 102, 51,255, 51, 0,204, 51, 0,255};
1368 
1369  TColor::CreateColorsCircle(kMagenta,"kMagenta",magenta);
1370  TColor::CreateColorsCircle(kRed, "kRed", red);
1371  TColor::CreateColorsCircle(kYellow, "kYellow", yellow);
1372  TColor::CreateColorsCircle(kGreen, "kGreen", green);
1373  TColor::CreateColorsCircle(kCyan, "kCyan", cyan);
1374  TColor::CreateColorsCircle(kBlue, "kBlue", blue);
1375 
1376  TColor::CreateColorsRectangle(kPink, "kPink", pink);
1377  TColor::CreateColorsRectangle(kOrange,"kOrange",orange);
1378  TColor::CreateColorsRectangle(kSpring,"kSpring",spring);
1379  TColor::CreateColorsRectangle(kTeal, "kTeal", teal);
1380  TColor::CreateColorsRectangle(kAzure, "kAzure", azure);
1381  TColor::CreateColorsRectangle(kViolet,"kViolet",violet);
1382 
1384 }
1385 
1386 ////////////////////////////////////////////////////////////////////////////////
1387 /// Static function returning the color number i in current palette.
1388 
1390 {
1391  Int_t ncolors = fgPalette.fN;
1392  if (ncolors == 0) return 0;
1393  Int_t icol = i%ncolors;
1394  if (icol < 0) icol = 0;
1395  return fgPalette.fArray[icol];
1396 }
1397 
1398 ////////////////////////////////////////////////////////////////////////////////
1399 /// Static function returning the current active palette.
1400 
1402 {
1403  return fgPalette;
1404 }
1405 
1406 ////////////////////////////////////////////////////////////////////////////////
1407 /// Static function returning number of colors in the color palette.
1408 
1410 {
1411  return fgPalette.fN;
1412 }
1413 
1414 ////////////////////////////////////////////////////////////////////////////////
1415 /// Static function returning kTRUE if some new colors have been defined after
1416 /// initialisation or since the last call to this method. This allows to avoid
1417 /// the colors and palette streaming in TCanvas::Streamer if not needed.
1418 
1420 {
1421  // After initialization gDefinedColors == 649. If it is bigger it means some new
1422  // colors have been defined
1423  Bool_t hasChanged = (gDefinedColors - gLastDefinedColors) > 50;
1425  return hasChanged;
1426 }
1427 
1428 ////////////////////////////////////////////////////////////////////////////////
1429 /// Return pixel value corresponding to this color. This pixel value can
1430 /// be used in the GUI classes. This call does not work in batch mode since
1431 /// it needs to communicate with the graphics system.
1432 
1434 {
1435  if (gVirtualX && !gROOT->IsBatch()) {
1436  if (gApplication) {
1439  }
1440  return gVirtualX->GetPixel(fNumber);
1441  }
1442 
1443  return 0;
1444 }
1445 
1446 ////////////////////////////////////////////////////////////////////////////////
1447 /// Static method to compute RGB from HLS. The l and s are between [0,1]
1448 /// and h is between [0,360]. The returned r,g,b triplet is between [0,1].
1449 
1450 void TColor::HLS2RGB(Float_t hue, Float_t light, Float_t satur,
1451  Float_t &r, Float_t &g, Float_t &b)
1452 {
1453 
1454  Float_t rh, rl, rs, rm1, rm2;
1455  rh = rl = rs = 0;
1456  if (hue > 0) { rh = hue; if (rh > 360) rh = 360; }
1457  if (light > 0) { rl = light; if (rl > 1) rl = 1; }
1458  if (satur > 0) { rs = satur; if (rs > 1) rs = 1; }
1459 
1460  if (rl <= 0.5)
1461  rm2 = rl*(1.0 + rs);
1462  else
1463  rm2 = rl + rs - rl*rs;
1464  rm1 = 2.0*rl - rm2;
1465 
1466  if (!rs) { r = rl; g = rl; b = rl; return; }
1467  r = HLStoRGB1(rm1, rm2, rh+120);
1468  g = HLStoRGB1(rm1, rm2, rh);
1469  b = HLStoRGB1(rm1, rm2, rh-120);
1470 }
1471 
1472 ////////////////////////////////////////////////////////////////////////////////
1473 /// Static method. Auxiliary to HLS2RGB().
1474 
1476 {
1477  Float_t hue = huei;
1478  if (hue > 360) hue = hue - 360;
1479  if (hue < 0) hue = hue + 360;
1480  if (hue < 60 ) return rn1 + (rn2-rn1)*hue/60;
1481  if (hue < 180) return rn2;
1482  if (hue < 240) return rn1 + (rn2-rn1)*(240-hue)/60;
1483  return rn1;
1484 }
1485 
1486 ////////////////////////////////////////////////////////////////////////////////
1487 /// Static method to compute RGB from HLS. The h,l,s are between [0,255].
1488 /// The returned r,g,b triplet is between [0,255].
1489 
1491 {
1492  Float_t hh, ll, ss, rr, gg, bb;
1493 
1494  hh = Float_t(h) * 360 / 255;
1495  ll = Float_t(l) / 255;
1496  ss = Float_t(s) / 255;
1497 
1498  TColor::HLStoRGB(hh, ll, ss, rr, gg, bb);
1499 
1500  r = (Int_t) (rr * 255);
1501  g = (Int_t) (gg * 255);
1502  b = (Int_t) (bb * 255);
1503 }
1504 
1505 ////////////////////////////////////////////////////////////////////////////////
1506 /// Static method to compute RGB from HSV:
1507 ///
1508 /// - The hue value runs from 0 to 360.
1509 /// - The saturation is the degree of strength or purity and is from 0 to 1.
1510 /// Purity is how much white is added to the color, so S=1 makes the purest
1511 /// color (no white).
1512 /// - Brightness value also ranges from 0 to 1, where 0 is the black.
1513 ///
1514 /// The returned r,g,b triplet is between [0,1].
1515 
1516 void TColor::HSV2RGB(Float_t hue, Float_t satur, Float_t value,
1517  Float_t &r, Float_t &g, Float_t &b)
1518 {
1519  Int_t i;
1520  Float_t f, p, q, t;
1521 
1522  if (satur==0) {
1523  // Achromatic (grey)
1524  r = g = b = value;
1525  return;
1526  }
1527 
1528  hue /= 60; // sector 0 to 5
1529  i = (Int_t)floor(hue);
1530  f = hue-i; // factorial part of hue
1531  p = value*(1-satur);
1532  q = value*(1-satur*f );
1533  t = value*(1-satur*(1-f));
1534 
1535  switch (i) {
1536  case 0:
1537  r = value;
1538  g = t;
1539  b = p;
1540  break;
1541  case 1:
1542  r = q;
1543  g = value;
1544  b = p;
1545  break;
1546  case 2:
1547  r = p;
1548  g = value;
1549  b = t;
1550  break;
1551  case 3:
1552  r = p;
1553  g = q;
1554  b = value;
1555  break;
1556  case 4:
1557  r = t;
1558  g = p;
1559  b = value;
1560  break;
1561  default:
1562  r = value;
1563  g = p;
1564  b = q;
1565  break;
1566  }
1567 }
1568 
1569 ////////////////////////////////////////////////////////////////////////////////
1570 /// List this color with its attributes.
1571 
1572 void TColor::ls(Option_t *) const
1573 {
1574  printf("Color:%d Red=%f Green=%f Blue=%f Alpha=%f Name=%s\n",
1576 }
1577 
1578 ////////////////////////////////////////////////////////////////////////////////
1579 /// Dump this color with its attributes.
1580 
1582 {
1583  ls();
1584 }
1585 
1586 ////////////////////////////////////////////////////////////////////////////////
1587 /// Static method to compute HLS from RGB. The r,g,b triplet is between
1588 /// [0,1], hue is between [0,360], light and satur are [0,1].
1589 
1591  Float_t &hue, Float_t &light, Float_t &satur)
1592 {
1593  Float_t rnorm, gnorm, bnorm, minval, maxval, msum, mdiff, r, g, b;
1594  minval = maxval =0 ;
1595  r = g = b = 0;
1596  if (rr > 0) { r = rr; if (r > 1) r = 1; }
1597  if (gg > 0) { g = gg; if (g > 1) g = 1; }
1598  if (bb > 0) { b = bb; if (b > 1) b = 1; }
1599 
1600  minval = r;
1601  if (g < minval) minval = g;
1602  if (b < minval) minval = b;
1603  maxval = r;
1604  if (g > maxval) maxval = g;
1605  if (b > maxval) maxval = b;
1606 
1607  rnorm = gnorm = bnorm = 0;
1608  mdiff = maxval - minval;
1609  msum = maxval + minval;
1610  light = 0.5 * msum;
1611  if (maxval != minval) {
1612  rnorm = (maxval - r)/mdiff;
1613  gnorm = (maxval - g)/mdiff;
1614  bnorm = (maxval - b)/mdiff;
1615  } else {
1616  satur = hue = 0;
1617  return;
1618  }
1619 
1620  if (light < 0.5)
1621  satur = mdiff/msum;
1622  else
1623  satur = mdiff/(2.0 - msum);
1624 
1625  if (r == maxval)
1626  hue = 60.0 * (6.0 + bnorm - gnorm);
1627  else if (g == maxval)
1628  hue = 60.0 * (2.0 + rnorm - bnorm);
1629  else
1630  hue = 60.0 * (4.0 + gnorm - rnorm);
1631 
1632  if (hue > 360)
1633  hue = hue - 360;
1634 }
1635 
1636 ////////////////////////////////////////////////////////////////////////////////
1637 /// Static method to compute HSV from RGB.
1638 ///
1639 /// - The input values:
1640 /// - r,g,b triplet is between [0,1].
1641 /// - The returned values:
1642 /// - The hue value runs from 0 to 360.
1643 /// - The saturation is the degree of strength or purity and is from 0 to 1.
1644 /// Purity is how much white is added to the color, so S=1 makes the purest
1645 /// color (no white).
1646 /// - Brightness value also ranges from 0 to 1, where 0 is the black.
1647 
1649  Float_t &hue, Float_t &satur, Float_t &value)
1650 {
1651  Float_t min, max, delta;
1652 
1653  min = TMath::Min(TMath::Min(r, g), b);
1654  max = TMath::Max(TMath::Max(r, g), b);
1655  value = max;
1656 
1657  delta = max - min;
1658 
1659  if (max != 0) {
1660  satur = delta/max;
1661  } else {
1662  satur = 0;
1663  hue = -1;
1664  return;
1665  }
1666 
1667  if (r == max) {
1668  hue = (g-b)/delta;
1669  } else if (g == max) {
1670  hue = 2+(b-r)/delta;
1671  } else {
1672  hue = 4+(r-g)/delta;
1673  }
1674 
1675  hue *= 60;
1676  if (hue < 0) hue += 360;
1677 }
1678 
1679 ////////////////////////////////////////////////////////////////////////////////
1680 /// Static method to compute HLS from RGB. The r,g,b triplet is between
1681 /// [0,255], hue, light and satur are between [0,255].
1682 
1684 {
1685  Float_t rr, gg, bb, hue, light, satur;
1686 
1687  rr = Float_t(r) / 255;
1688  gg = Float_t(g) / 255;
1689  bb = Float_t(b) / 255;
1690 
1691  TColor::RGBtoHLS(rr, gg, bb, hue, light, satur);
1692 
1693  h = (Int_t) (hue/360 * 255);
1694  l = (Int_t) (light * 255);
1695  s = (Int_t) (satur * 255);
1696 }
1697 
1698 ////////////////////////////////////////////////////////////////////////////////
1699 /// Initialize this color and its associated colors.
1700 
1702 {
1704  fRed = r;
1705  fGreen = g;
1706  fBlue = b;
1707 
1708  if (fRed < 0) return;
1709 
1711 
1712  Int_t nplanes = 16;
1713  if (gVirtualX) gVirtualX->GetPlanes(nplanes);
1714  if (nplanes == 0) nplanes = 16;
1715 
1716  // allocate color now (can be delayed when we have a large colormap)
1717 #ifndef R__WIN32
1718  if (nplanes < 15)
1719 #endif
1720  Allocate();
1721 
1722  if (fNumber > 50) return;
1723 
1724  // now define associated colors for WBOX shading
1725  Float_t dr, dg, db, lr, lg, lb;
1726 
1727  // set dark color
1728  HLStoRGB(fHue, 0.7*fLight, fSaturation, dr, dg, db);
1729  TColor *dark = gROOT->GetColor(100+fNumber);
1730  if (dark) {
1731  if (nplanes > 8) dark->SetRGB(dr, dg, db);
1732  else dark->SetRGB(0.3,0.3,0.3);
1733  }
1734 
1735  // set light color
1736  HLStoRGB(fHue, 1.2*fLight, fSaturation, lr, lg, lb);
1737  TColor *light = gROOT->GetColor(150+fNumber);
1738  if (light) {
1739  if (nplanes > 8) light->SetRGB(lr, lg, lb);
1740  else light->SetRGB(0.8,0.8,0.8);
1741  }
1742  gDefinedColors++;
1743 }
1744 
1745 ////////////////////////////////////////////////////////////////////////////////
1746 /// Make this color known to the graphics system.
1747 
1749 {
1750  if (gVirtualX && !gROOT->IsBatch())
1751 
1752  gVirtualX->SetRGB(fNumber, GetRed(), GetGreen(), GetBlue());
1753 }
1754 
1755 ////////////////////////////////////////////////////////////////////////////////
1756 /// Static method returning color number for color specified by
1757 /// hex color string of form: "#rrggbb", where rr, gg and bb are in
1758 /// hex between [0,FF], e.g. "#c0c0c0".
1759 ///
1760 /// The color retrieval is done using a threshold defined by SetColorThreshold.
1761 ///
1762 /// If specified color does not exist it will be created with as
1763 /// name "#rrggbb" with rr, gg and bb in hex between [0,FF].
1764 
1765 Int_t TColor::GetColor(const char *hexcolor)
1766 {
1767  if (hexcolor && *hexcolor == '#') {
1768  Int_t r, g, b;
1769  if (sscanf(hexcolor+1, "%02x%02x%02x", &r, &g, &b) == 3)
1770  return GetColor(r, g, b);
1771  }
1772  ::Error("TColor::GetColor(const char*)", "incorrect color string");
1773  return 0;
1774 }
1775 
1776 ////////////////////////////////////////////////////////////////////////////////
1777 /// Static method returning color number for color specified by
1778 /// r, g and b. The r,g,b should be in the range [0,1].
1779 ///
1780 /// The color retrieval is done using a threshold defined by SetColorThreshold.
1781 ///
1782 /// If specified color does not exist it will be created
1783 /// with as name "#rrggbb" with rr, gg and bb in hex between
1784 /// [0,FF].
1785 
1787 {
1788  Int_t rr, gg, bb;
1789  rr = Int_t(r * 255);
1790  gg = Int_t(g * 255);
1791  bb = Int_t(b * 255);
1792 
1793  return GetColor(rr, gg, bb);
1794 }
1795 
1796 ////////////////////////////////////////////////////////////////////////////////
1797 /// Static method returning color number for color specified by
1798 /// system dependent pixel value. Pixel values can be obtained, e.g.,
1799 /// from the GUI color picker.
1800 ///
1801 /// The color retrieval is done using a threshold defined by SetColorThreshold.
1802 
1803 
1805 {
1806  Int_t r, g, b;
1807 
1808  Pixel2RGB(pixel, r, g, b);
1809 
1810  return GetColor(r, g, b);
1811 }
1812 
1813 ////////////////////////////////////////////////////////////////////////////////
1814 /// This method specifies the color threshold used by GetColor to retrieve a color.
1815 ///
1816 /// \param[in] t Color threshold. By default is equal to 1./31. or 1./255.
1817 /// depending on the number of available color planes.
1818 ///
1819 /// When GetColor is called, it scans the defined colors and compare them to the
1820 /// requested color.
1821 /// If the Red Green and Blue values passed to GetColor are Rr Gr Br
1822 /// and Rd Gd Bd the values of a defined color. These two colors are considered equal
1823 /// if (abs(Rr-Rd) < t & abs(Br-Bd) < t & abs(Br-Bd) < t). If this test passes,
1824 /// the color defined by Rd Gd Bd is returned by GetColor.
1825 ///
1826 /// To make sure GetColor will return a color having exactly the requested
1827 /// R G B values it is enough to specify a nul :
1828 /// ~~~ {.cpp}
1829 /// TColor::SetColorThreshold(0.);
1830 /// ~~~
1831 ///
1832 /// To reset the color threshold to its default value it is enough to do:
1833 /// ~~~ {.cpp}
1834 /// TColor::SetColorThreshold(-1.);
1835 /// ~~~
1836 
1838 {
1839  gColorThreshold = t;
1840 }
1841 
1842 ////////////////////////////////////////////////////////////////////////////////
1843 /// Static method returning color number for color specified by
1844 /// r, g and b. The r,g,b should be in the range [0,255].
1845 /// If the specified color does not exist it will be created
1846 /// with as name "#rrggbb" with rr, gg and bb in hex between
1847 /// [0,FF].
1848 ///
1849 /// The color retrieval is done using a threshold defined by SetColorThreshold.
1850 
1851 
1853 {
1855  if (r < 0) r = 0;
1856  if (g < 0) g = 0;
1857  if (b < 0) b = 0;
1858  if (r > 255) r = 255;
1859  if (g > 255) g = 255;
1860  if (b > 255) b = 255;
1861 
1862  // Get list of all defined colors
1863  TObjArray *colors = (TObjArray*) gROOT->GetListOfColors();
1864 
1865  TColor *color = 0;
1866 
1867  // Look for color by name
1868  if ((color = (TColor*) colors->FindObject(Form("#%02x%02x%02x", r, g, b))))
1869  // We found the color by name, so we use that right away
1870  return color->GetNumber();
1871 
1872  Float_t rr, gg, bb;
1873  rr = Float_t(r)/255.;
1874  gg = Float_t(g)/255.;
1875  bb = Float_t(b)/255.;
1876 
1877  TIter next(colors);
1878 
1879  Float_t thres;
1880  if (gColorThreshold >= 0) {
1881  thres = gColorThreshold;
1882  } else {
1883  Int_t nplanes = 16;
1884  thres = 1.0/31.0; // 5 bits per color : 0 - 0x1F !
1885  if (gVirtualX) gVirtualX->GetPlanes(nplanes);
1886  if (nplanes >= 24) thres = 1.0/255.0; // 8 bits per color : 0 - 0xFF !
1887  }
1888 
1889  // Loop over all defined colors
1890  while ((color = (TColor*)next())) {
1891  if (TMath::Abs(color->GetRed() - rr) > thres) continue;
1892  if (TMath::Abs(color->GetGreen() - gg) > thres) continue;
1893  if (TMath::Abs(color->GetBlue() - bb) > thres) continue;
1894  // We found a matching color in the color table
1895  return color->GetNumber();
1896  }
1897 
1898  // We didn't find a matching color in the color table, so we
1899  // add it. Note name is of the form "#rrggbb" where rr, etc. are
1900  // hexadecimal numbers.
1901  color = new TColor(colors->GetLast()+1, rr, gg, bb,
1902  Form("#%02x%02x%02x", r, g, b));
1903 
1904  return color->GetNumber();
1905 }
1906 
1907 ////////////////////////////////////////////////////////////////////////////////
1908 /// Static function: Returns the bright color number corresponding to n
1909 /// If the TColor object does not exist, it is created.
1910 /// The convention is that the bright color nb = n+150
1911 
1913 {
1914  if (n < 0) return -1;
1915 
1916  // Get list of all defined colors
1917  TObjArray *colors = (TObjArray*) gROOT->GetListOfColors();
1918  Int_t ncolors = colors->GetSize();
1919  // Get existing color at index n
1920  TColor *color = 0;
1921  if (n < ncolors) color = (TColor*)colors->At(n);
1922  if (!color) return -1;
1923 
1924  //Get the rgb of the the new bright color corresponding to color n
1925  Float_t r,g,b;
1926  HLStoRGB(color->GetHue(), 1.2*color->GetLight(), color->GetSaturation(), r, g, b);
1927 
1928  //Build the bright color (unless the slot nb is already used)
1929  Int_t nb = n+150;
1930  TColor *colorb = 0;
1931  if (nb < ncolors) colorb = (TColor*)colors->At(nb);
1932  if (colorb) return nb;
1933  colorb = new TColor(nb,r,g,b);
1934  colorb->SetName(Form("%s_bright",color->GetName()));
1935  colors->AddAtAndExpand(colorb,nb);
1936  return nb;
1937 }
1938 
1939 ////////////////////////////////////////////////////////////////////////////////
1940 /// Static function: Returns the dark color number corresponding to n
1941 /// If the TColor object does not exist, it is created.
1942 /// The convention is that the dark color nd = n+100
1943 
1945 {
1946  if (n < 0) return -1;
1947 
1948  // Get list of all defined colors
1949  TObjArray *colors = (TObjArray*) gROOT->GetListOfColors();
1950  Int_t ncolors = colors->GetSize();
1951  // Get existing color at index n
1952  TColor *color = 0;
1953  if (n < ncolors) color = (TColor*)colors->At(n);
1954  if (!color) return -1;
1955 
1956  //Get the rgb of the the new dark color corresponding to color n
1957  Float_t r,g,b;
1958  HLStoRGB(color->GetHue(), 0.7*color->GetLight(), color->GetSaturation(), r, g, b);
1959 
1960  //Build the dark color (unless the slot nd is already used)
1961  Int_t nd = n+100;
1962  TColor *colord = 0;
1963  if (nd < ncolors) colord = (TColor*)colors->At(nd);
1964  if (colord) return nd;
1965  colord = new TColor(nd,r,g,b);
1966  colord->SetName(Form("%s_dark",color->GetName()));
1967  colors->AddAtAndExpand(colord,nd);
1968  return nd;
1969 }
1970 
1971 ////////////////////////////////////////////////////////////////////////////////
1972 /// Static function: Returns the transparent color number corresponding to n.
1973 /// The transparency level is given by the alpha value a.
1974 
1976 {
1977  if (n < 0) return -1;
1978 
1979  TColor *color = gROOT->GetColor(n);
1980  if (color) {
1981  TColor *colort = new TColor(gROOT->GetListOfColors()->GetLast()+1,
1982  color->GetRed(), color->GetGreen(), color->GetBlue());
1983  colort->SetAlpha(a);
1984  colort->SetName(Form("%s_transparent",color->GetName()));
1985  return colort->GetNumber();
1986  } else {
1987  ::Error("TColor::GetColorTransparent", "color with index %d not defined", n);
1988  return -1;
1989  }
1990 }
1991 
1992 ////////////////////////////////////////////////////////////////////////////////
1993 /// Static function: Returns a free color index which can be used to define
1994 /// a user custom color.
1995 ///
1996 /// ~~~ {.cpp}
1997 /// Int_t ci = TColor::GetFreeColorIndex();
1998 /// TColor *color = new TColor(ci, 0.1, 0.2, 0.3);
1999 /// ~~~
2000 
2002 {
2003  return gHighestColorIndex+1;
2004 }
2005 
2006 ////////////////////////////////////////////////////////////////////////////////
2007 /// Static method that given a color index number, returns the corresponding
2008 /// pixel value. This pixel value can be used in the GUI classes. This call
2009 /// does not work in batch mode since it needs to communicate with the
2010 /// graphics system.
2011 
2013 {
2015  TColor *color = gROOT->GetColor(ci);
2016  if (color)
2017  return color->GetPixel();
2018  else
2019  ::Warning("TColor::Number2Pixel", "color with index %d not defined", ci);
2020 
2021  return 0;
2022 }
2023 
2024 ////////////////////////////////////////////////////////////////////////////////
2025 /// Convert r,g,b to graphics system dependent pixel value.
2026 /// The r,g,b triplet must be [0,1].
2027 
2029 {
2030  if (r < 0) r = 0;
2031  if (g < 0) g = 0;
2032  if (b < 0) b = 0;
2033  if (r > 1) r = 1;
2034  if (g > 1) g = 1;
2035  if (b > 1) b = 1;
2036 
2037  ColorStruct_t color;
2038  color.fRed = UShort_t(r * 65535);
2039  color.fGreen = UShort_t(g * 65535);
2040  color.fBlue = UShort_t(b * 65535);
2041  color.fMask = kDoRed | kDoGreen | kDoBlue;
2042  gVirtualX->AllocColor(gVirtualX->GetColormap(), color);
2043  return color.fPixel;
2044 }
2045 
2046 ////////////////////////////////////////////////////////////////////////////////
2047 /// Convert r,g,b to graphics system dependent pixel value.
2048 /// The r,g,b triplet must be [0,255].
2049 
2051 {
2052  if (r < 0) r = 0;
2053  if (g < 0) g = 0;
2054  if (b < 0) b = 0;
2055  if (r > 255) r = 255;
2056  if (g > 255) g = 255;
2057  if (b > 255) b = 255;
2058 
2059  ColorStruct_t color;
2060  color.fRed = UShort_t(r * 257); // 65535/255
2061  color.fGreen = UShort_t(g * 257);
2062  color.fBlue = UShort_t(b * 257);
2063  color.fMask = kDoRed | kDoGreen | kDoBlue;
2064  gVirtualX->AllocColor(gVirtualX->GetColormap(), color);
2065  return color.fPixel;
2066 }
2067 
2068 ////////////////////////////////////////////////////////////////////////////////
2069 /// Convert machine dependent pixel value (obtained via RGB2Pixel or
2070 /// via Number2Pixel() or via TColor::GetPixel()) to r,g,b triplet.
2071 /// The r,g,b triplet will be [0,1].
2072 
2074 {
2075  ColorStruct_t color;
2076  color.fPixel = pixel;
2077  gVirtualX->QueryColor(gVirtualX->GetColormap(), color);
2078  r = (Float_t)color.fRed / 65535;
2079  g = (Float_t)color.fGreen / 65535;
2080  b = (Float_t)color.fBlue / 65535;
2081 }
2082 
2083 ////////////////////////////////////////////////////////////////////////////////
2084 /// Convert machine dependent pixel value (obtained via RGB2Pixel or
2085 /// via Number2Pixel() or via TColor::GetPixel()) to r,g,b triplet.
2086 /// The r,g,b triplet will be [0,255].
2087 
2089 {
2090  ColorStruct_t color;
2091  color.fPixel = pixel;
2092  gVirtualX->QueryColor(gVirtualX->GetColormap(), color);
2093  r = color.fRed / 257;
2094  g = color.fGreen / 257;
2095  b = color.fBlue / 257;
2096 }
2097 
2098 ////////////////////////////////////////////////////////////////////////////////
2099 /// Convert machine dependent pixel value (obtained via RGB2Pixel or
2100 /// via Number2Pixel() or via TColor::GetPixel()) to a hexadecimal string.
2101 /// This string can be directly passed to, for example,
2102 /// TGClient::GetColorByName(). String will be reused so copy immediately
2103 /// if needed.
2104 
2106 {
2107  static TString tempbuf;
2108  Int_t r, g, b;
2109  Pixel2RGB(pixel, r, g, b);
2110  tempbuf.Form("#%02x%02x%02x", r, g, b);
2111  return tempbuf;
2112 }
2113 
2114 ////////////////////////////////////////////////////////////////////////////////
2115 /// Save a color with index > 228 as a C++ statement(s) on output stream out.
2116 
2117 void TColor::SaveColor(std::ostream &out, Int_t ci)
2118 {
2119  char quote = '"';
2120  Float_t r,g,b,a;
2121  Int_t ri, gi, bi;
2122  TString cname;
2123 
2124  TColor *c = gROOT->GetColor(ci);
2125  if (c) {
2126  c->GetRGB(r, g, b);
2127  a = c->GetAlpha();
2128  } else {
2129  return;
2130  }
2131 
2132  if (gROOT->ClassSaved(TColor::Class())) {
2133  out << std::endl;
2134  } else {
2135  out << std::endl;
2136  out << " Int_t ci; // for color index setting" << std::endl;
2137  out << " TColor *color; // for color definition with alpha" << std::endl;
2138  }
2139 
2140  if (a<1) {
2141  out<<" ci = "<<ci<<";"<<std::endl;
2142  out<<" color = new TColor(ci, "<<r<<", "<<g<<", "<<b<<", "
2143  <<"\" \", "<<a<<");"<<std::endl;
2144  } else {
2145  ri = (Int_t)(255*r);
2146  gi = (Int_t)(255*g);
2147  bi = (Int_t)(255*b);
2148  cname.Form("#%02x%02x%02x", ri, gi, bi);
2149  out<<" ci = TColor::GetColor("<<quote<<cname.Data()<<quote<<");"<<std::endl;
2150  }
2151 }
2152 
2153 ////////////////////////////////////////////////////////////////////////////////
2154 /// Return whether all colors return grayscale values.
2156 {
2157  return fgGrayscaleMode;
2158 }
2159 
2160 ////////////////////////////////////////////////////////////////////////////////
2161 /// Set whether all colors should return grayscale values.
2162 
2163 void TColor::SetGrayscale(Bool_t set /*= kTRUE*/)
2164 {
2165  if (fgGrayscaleMode == set) return;
2166 
2167  fgGrayscaleMode = set;
2168 
2169  if (!gVirtualX || gROOT->IsBatch()) return;
2170 
2172  TIter iColor(gROOT->GetListOfColors());
2173  TColor* color = 0;
2174  while ((color = (TColor*) iColor()))
2175  color->Allocate();
2176 }
2177 
2178 ////////////////////////////////////////////////////////////////////////////////
2179 /// Static function creating a color table with several connected linear gradients.
2180 ///
2181 /// - Number: The number of end point colors that will form the gradients.
2182 /// Must be at least 2.
2183 /// - Stops: Where in the whole table the end point colors should lie.
2184 /// Each entry must be on [0, 1], each entry must be greater than
2185 /// the previous entry.
2186 /// - Red, Green, Blue: The end point color values.
2187 /// Each entry must be on [0, 1]
2188 /// - NColors: Total number of colors in the table. Must be at least 1.
2189 ///
2190 /// Returns a positive value on success and -1 on error.
2191 ///
2192 /// The table is constructed by tracing lines between the given points in
2193 /// RGB space. Each color value may have a value between 0 and 1. The
2194 /// difference between consecutive "Stops" values gives the fraction of
2195 /// space in the whole table that should be used for the interval between
2196 /// the corresponding color values.
2197 ///
2198 /// Normally the first element of Stops should be 0 and the last should be 1.
2199 /// If this is not true, fewer than NColors will be used in proportion with
2200 /// the total interval between the first and last elements of Stops.
2201 ///
2202 /// This definition is similar to the povray-definition of gradient
2203 /// color tables.
2204 ///
2205 /// For instance:
2206 /// ~~~ {.cpp}
2207 /// UInt_t Number = 3;
2208 /// Double_t Red[3] = { 0.0, 1.0, 1.0 };
2209 /// Double_t Green[3] = { 0.0, 0.0, 1.0 };
2210 /// Double_t Blue[3] = { 1.0, 0.0, 1.0 };
2211 /// Double_t Stops[3] = { 0.0, 0.4, 1.0 };
2212 /// ~~~
2213 /// This defines a table in which there are three color end points:
2214 /// RGB = {0, 0, 1}, {1, 0, 0}, and {1, 1, 1} = blue, red, white
2215 /// The first 40% of the table is used to go linearly from blue to red.
2216 /// The remaining 60% of the table is used to go linearly from red to white.
2217 ///
2218 /// If you define a very short interval such that less than one color fits
2219 /// in it, no colors at all will be allocated. If this occurs for all
2220 /// intervals, ROOT will revert to the default palette.
2221 ///
2222 /// Original code by Andreas Zoglauer (zog@mpe.mpg.de)
2223 
2225  Double_t* Red, Double_t* Green,
2226  Double_t* Blue, UInt_t NColors, Float_t alpha)
2227 {
2229 
2230  UInt_t g, c;
2231  UInt_t nPalette = 0;
2232  Int_t *palette = new Int_t[NColors+1];
2233  UInt_t nColorsGradient;
2234 
2235  if(Number < 2 || NColors < 1){
2236  delete [] palette;
2237  return -1;
2238  }
2239 
2240  // Check if all RGB values are between 0.0 and 1.0 and
2241  // Stops goes from 0.0 to 1.0 in increasing order.
2242  for (c = 0; c < Number; c++) {
2243  if (Red[c] < 0 || Red[c] > 1.0 ||
2244  Green[c] < 0 || Green[c] > 1.0 ||
2245  Blue[c] < 0 || Blue[c] > 1.0 ||
2246  Stops[c] < 0 || Stops[c] > 1.0) {
2247  delete [] palette;
2248  return -1;
2249  }
2250  if (c >= 1) {
2251  if (Stops[c-1] > Stops[c]) {
2252  delete [] palette;
2253  return -1;
2254  }
2255  }
2256  }
2257 
2258  // Now create the colors and add them to the default palette:
2259 
2260  // For each defined gradient...
2261  for (g = 1; g < Number; g++) {
2262  // create the colors...
2263  nColorsGradient = (Int_t) (floor(NColors*Stops[g]) - floor(NColors*Stops[g-1]));
2264  for (c = 0; c < nColorsGradient; c++) {
2265  new TColor( Float_t(Red[g-1] + c * (Red[g] - Red[g-1]) / nColorsGradient),
2266  Float_t(Green[g-1] + c * (Green[g] - Green[g-1])/ nColorsGradient),
2267  Float_t(Blue[g-1] + c * (Blue[g] - Blue[g-1]) / nColorsGradient),
2268  alpha);
2269  palette[nPalette] = gHighestColorIndex;
2270  nPalette++;
2271  }
2272  }
2273 
2274  TColor::SetPalette(nPalette, palette);
2275  delete [] palette;
2276  return gHighestColorIndex + 1 - NColors;
2277 }
2278 
2279 
2280 ////////////////////////////////////////////////////////////////////////////////
2281 /// Static function.
2282 /// The color palette is used by the histogram classes
2283 /// (see TH1::Draw options).
2284 /// For example TH1::Draw("col") draws a 2-D histogram with cells
2285 /// represented by a box filled with a color CI function of the cell content.
2286 /// if the cell content is N, the color CI used will be the color number
2287 /// in colors[N],etc. If the maximum cell content is > ncolors, all
2288 /// cell contents are scaled to ncolors.
2289 ///
2290 /// `if ncolors <= 0` a default palette (see below) of 50 colors is
2291 /// defined. The colors defined in this palette are OK for coloring pads, labels.
2292 ///
2293 /// ~~~ {.cpp}
2294 /// index 0->9 : grey colors from light to dark grey
2295 /// index 10->19 : "brown" colors
2296 /// index 20->29 : "blueish" colors
2297 /// index 30->39 : "redish" colors
2298 /// index 40->49 : basic colors
2299 /// ~~~
2300 ///
2301 /// `if ncolors == 1 && colors == 0`, a Rainbow Color map is created
2302 /// with 50 colors. It is kept for backward compatibility. Better palettes like
2303 /// kBird are recommended.
2304 ///
2305 /// High quality predefined palettes with 255 colors are available when `colors == 0`.
2306 /// The following value of `ncolors` give access to:
2307 ///
2308 /// ~~~ {.cpp}
2309 /// if ncolors = 51 and colors=0, a Deep Sea palette is used.
2310 /// if ncolors = 52 and colors=0, a Grey Scale palette is used.
2311 /// if ncolors = 53 and colors=0, a Dark Body Radiator palette is used.
2312 /// if ncolors = 54 and colors=0, a Two-Color Hue palette is used.(dark blue through neutral gray to bright yellow)
2313 /// if ncolors = 55 and colors=0, a Rain Bow palette is used.
2314 /// if ncolors = 56 and colors=0, an Inverted Dark Body Radiator palette is used.
2315 /// if ncolors = 57 and colors=0, a monotonically increasing L value palette is used.
2316 /// if ncolors = 58 and colors=0, a Cubehelix palette is used
2317 /// (Cf. Dave Green's "cubehelix" colour scheme at http://www.mrao.cam.ac.uk/~dag/CUBEHELIX/)
2318 /// if ncolors = 59 and colors=0, a Green Red Violet palette is used.
2319 /// if ncolors = 60 and colors=0, a Blue Red Yellow palette is used.
2320 /// if ncolors = 61 and colors=0, an Ocean palette is used.
2321 /// if ncolors = 62 and colors=0, a Color Printable On Grey palette is used.
2322 /// if ncolors = 63 and colors=0, an Alpine palette is used.
2323 /// if ncolors = 64 and colors=0, an Aquamarine palette is used.
2324 /// if ncolors = 65 and colors=0, an Army palette is used.
2325 /// if ncolors = 66 and colors=0, an Atlantic palette is used.
2326 /// if ncolors = 67 and colors=0, an Aurora palette is used.
2327 /// if ncolors = 68 and colors=0, an Avocado palette is used.
2328 /// if ncolors = 69 and colors=0, a Beach palette is used.
2329 /// if ncolors = 70 and colors=0, a Black Body palette is used.
2330 /// if ncolors = 71 and colors=0, a Blue Green Yellow palette is used.
2331 /// if ncolors = 72 and colors=0, a Brown Cyan palette is used.
2332 /// if ncolors = 73 and colors=0, a CMYK palette is used.
2333 /// if ncolors = 74 and colors=0, a Candy palette is used.
2334 /// if ncolors = 75 and colors=0, a Cherry palette is used.
2335 /// if ncolors = 76 and colors=0, a Coffee palette is used.
2336 /// if ncolors = 77 and colors=0, a Dark Rain Bow palette is used.
2337 /// if ncolors = 78 and colors=0, a Dark Terrain palette is used.
2338 /// if ncolors = 79 and colors=0, a Fall palette is used.
2339 /// if ncolors = 80 and colors=0, a Fruit Punch palette is used.
2340 /// if ncolors = 81 and colors=0, a Fuchsia palette is used.
2341 /// if ncolors = 82 and colors=0, a Grey Yellow palette is used.
2342 /// if ncolors = 83 and colors=0, a Green Brown Terrain palette is used.
2343 /// if ncolors = 84 and colors=0, a Green Pink palette is used.
2344 /// if ncolors = 85 and colors=0, an Island palette is used.
2345 /// if ncolors = 86 and colors=0, a Lake palette is used.
2346 /// if ncolors = 87 and colors=0, a Light Temperature palette is used.
2347 /// if ncolors = 88 and colors=0, a Light Terrain palette is used.
2348 /// if ncolors = 89 and colors=0, a Mint palette is used.
2349 /// if ncolors = 90 and colors=0, a Neon palette is used.
2350 /// if ncolors = 91 and colors=0, a Pastel palette is used.
2351 /// if ncolors = 92 and colors=0, a Pearl palette is used.
2352 /// if ncolors = 93 and colors=0, a Pigeon palette is used.
2353 /// if ncolors = 94 and colors=0, a Plum palette is used.
2354 /// if ncolors = 95 and colors=0, a Red Blue palette is used.
2355 /// if ncolors = 96 and colors=0, a Rose palette is used.
2356 /// if ncolors = 97 and colors=0, a Rust palette is used.
2357 /// if ncolors = 98 and colors=0, a Sandy Terrain palette is used.
2358 /// if ncolors = 99 and colors=0, a Sienna palette is used.
2359 /// if ncolors = 100 and colors=0, a Solar palette is used.
2360 /// if ncolors = 101 and colors=0, a South West palette is used.
2361 /// if ncolors = 102 and colors=0, a Starry Night palette is used.
2362 /// if ncolors = 103 and colors=0, a Sunset palette is used.
2363 /// if ncolors = 104 and colors=0, a Temperature Map palette is used.
2364 /// if ncolors = 105 and colors=0, a Thermometer palette is used.
2365 /// if ncolors = 106 and colors=0, a Valentine palette is used.
2366 /// if ncolors = 107 and colors=0, a Visible Spectrum palette is used.
2367 /// if ncolors = 108 and colors=0, a Water Melon palette is used.
2368 /// if ncolors = 109 and colors=0, a Cool palette is used.
2369 /// if ncolors = 110 and colors=0, a Copper palette is used.
2370 /// if ncolors = 111 and colors=0, a Gist Earth palette is used.
2371 /// if ncolors = 112 and colors=0, a Viridis palette is used.
2372 /// if ncolors = 113 and colors=0, a Cividis palette is used.
2373 /// ~~~
2374 /// These palettes can also be accessed by names:
2375 /// ~~~ {.cpp}
2376 /// kDeepSea=51, kGreyScale=52, kDarkBodyRadiator=53,
2377 /// kBlueYellow= 54, kRainBow=55, kInvertedDarkBodyRadiator=56,
2378 /// kBird=57, kCubehelix=58, kGreenRedViolet=59,
2379 /// kBlueRedYellow=60, kOcean=61, kColorPrintableOnGrey=62,
2380 /// kAlpine=63, kAquamarine=64, kArmy=65,
2381 /// kAtlantic=66, kAurora=67, kAvocado=68,
2382 /// kBeach=69, kBlackBody=70, kBlueGreenYellow=71,
2383 /// kBrownCyan=72, kCMYK=73, kCandy=74,
2384 /// kCherry=75, kCoffee=76, kDarkRainBow=77,
2385 /// kDarkTerrain=78, kFall=79, kFruitPunch=80,
2386 /// kFuchsia=81, kGreyYellow=82, kGreenBrownTerrain=83,
2387 /// kGreenPink=84, kIsland=85, kLake=86,
2388 /// kLightTemperature=87, kLightTerrain=88, kMint=89,
2389 /// kNeon=90, kPastel=91, kPearl=92,
2390 /// kPigeon=93, kPlum=94, kRedBlue=95,
2391 /// kRose=96, kRust=97, kSandyTerrain=98,
2392 /// kSienna=99, kSolar=100, kSouthWest=101,
2393 /// kStarryNight=102, kSunset=103, kTemperatureMap=104,
2394 /// kThermometer=105, kValentine=106, kVisibleSpectrum=107,
2395 /// kWaterMelon=108, kCool=109, kCopper=110,
2396 /// kGistEarth=111 kViridis=112, kCividis=113
2397 /// ~~~
2398 /// For example:
2399 /// ~~~ {.cpp}
2400 /// gStyle->SetPalette(kBird);
2401 /// ~~~
2402 /// Set the current palette as "Bird" (number 57).
2403 ///
2404 /// The color numbers specified in the palette can be viewed by selecting
2405 /// the item "colors" in the "VIEW" menu of the canvas toolbar.
2406 /// The color parameters can be changed via TColor::SetRGB.
2407 ///
2408 /// Note that when drawing a 2D histogram `h2` with the option "COL" or
2409 /// "COLZ" or with any "CONT" options using the color map, the number of colors
2410 /// used is defined by the number of contours `n` specified with:
2411 /// `h2->SetContour(n)`
2412 
2414 {
2415  Int_t i;
2416 
2417  static Int_t paletteType = 0;
2418 
2419  Int_t palette[50] = {19,18,17,16,15,14,13,12,11,20,
2420  21,22,23,24,25,26,27,28,29,30, 8,
2421  31,32,33,34,35,36,37,38,39,40, 9,
2422  41,42,43,44,45,47,48,49,46,50, 2,
2423  7, 6, 5, 4, 3, 2,1};
2424 
2425  // set default palette (pad type)
2426  if (ncolors <= 0) {
2427  ncolors = 50;
2428  fgPalette.Set(ncolors);
2429  for (i=0;i<ncolors;i++) fgPalette.fArray[i] = palette[i];
2430  paletteType = 1;
2431  return;
2432  }
2433 
2434  // set Rainbow Color map. Kept for backward compatibility.
2435  if (ncolors == 1 && colors == 0) {
2436  ncolors = 50;
2437  fgPalette.Set(ncolors);
2438  for (i=0;i<ncolors-1;i++) fgPalette.fArray[i] = 51+i;
2439  fgPalette.fArray[ncolors-1] = kRed; // the last color of this palette is red
2440  paletteType = 2;
2441  return;
2442  }
2443 
2444  // High quality palettes (255 levels)
2445  if (colors == 0 && ncolors>50) {
2446 
2447  if (!fgPalettesList.fN) fgPalettesList.Set(63); // Right now 63 high quality palettes
2448  Int_t Idx = (Int_t)fgPalettesList.fArray[ncolors-51]; // High quality palettes indices start at 51
2449 
2450  // This high quality palette has already been created. Reuse it.
2451  if (Idx > 0) {
2452  Double_t alphas = 10*(fgPalettesList.fArray[ncolors-51]-Idx);
2453  Bool_t same_alpha = TMath::Abs(alpha-alphas) < 0.0001;
2454  if (paletteType == ncolors && same_alpha) return; // The current palette is already this one.
2455  fgPalette.Set(255); // High quality palettes have 255 entries
2456  for (i=0;i<255;i++) fgPalette.fArray[i] = Idx+i;
2457  paletteType = ncolors;
2458 
2459  // restore the palette transparency if needed
2460  if (alphas>0 && !same_alpha) {
2461  TColor *ca;
2462  for (i=0;i<255;i++) {
2463  ca = gROOT->GetColor(Idx+i);
2464  ca->SetAlpha(alpha);
2465  }
2466  fgPalettesList.fArray[paletteType-51] = (Double_t)Idx+alpha/10.;
2467  }
2468  return;
2469  }
2470 
2472  Double_t stops[9] = { 0.0000, 0.1250, 0.2500, 0.3750, 0.5000, 0.6250, 0.7500, 0.8750, 1.0000};
2473 
2474  switch (ncolors) {
2475  // Deep Sea
2476  case 51:
2477  {
2478  Double_t red[9] = { 0./255., 9./255., 13./255., 17./255., 24./255., 32./255., 27./255., 25./255., 29./255.};
2479  Double_t green[9] = { 0./255., 0./255., 0./255., 2./255., 37./255., 74./255., 113./255., 160./255., 221./255.};
2480  Double_t blue[9] = { 28./255., 42./255., 59./255., 78./255., 98./255., 129./255., 154./255., 184./255., 221./255.};
2481  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2482  }
2483  break;
2484 
2485  // Grey Scale
2486  case 52:
2487  {
2488  Double_t red[9] = { 0./255., 32./255., 64./255., 96./255., 128./255., 160./255., 192./255., 224./255., 255./255.};
2489  Double_t green[9] = { 0./255., 32./255., 64./255., 96./255., 128./255., 160./255., 192./255., 224./255., 255./255.};
2490  Double_t blue[9] = { 0./255., 32./255., 64./255., 96./255., 128./255., 160./255., 192./255., 224./255., 255./255.};
2491  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2492  }
2493  break;
2494 
2495  // Dark Body Radiator
2496  case 53:
2497  {
2498  Double_t red[9] = { 0./255., 45./255., 99./255., 156./255., 212./255., 230./255., 237./255., 234./255., 242./255.};
2499  Double_t green[9] = { 0./255., 0./255., 0./255., 45./255., 101./255., 168./255., 238./255., 238./255., 243./255.};
2500  Double_t blue[9] = { 0./255., 1./255., 1./255., 3./255., 9./255., 8./255., 11./255., 95./255., 230./255.};
2501  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2502  }
2503  break;
2504 
2505  // Two-color hue (dark blue through neutral gray to bright yellow)
2506  case 54:
2507  {
2508  Double_t red[9] = { 0./255., 22./255., 44./255., 68./255., 93./255., 124./255., 160./255., 192./255., 237./255.};
2509  Double_t green[9] = { 0./255., 16./255., 41./255., 67./255., 93./255., 125./255., 162./255., 194./255., 241./255.};
2510  Double_t blue[9] = { 97./255., 100./255., 99./255., 99./255., 93./255., 68./255., 44./255., 26./255., 74./255.};
2511  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2512  }
2513  break;
2514 
2515  // Rain Bow
2516  case 55:
2517  {
2518  Double_t red[9] = { 0./255., 5./255., 15./255., 35./255., 102./255., 196./255., 208./255., 199./255., 110./255.};
2519  Double_t green[9] = { 0./255., 48./255., 124./255., 192./255., 206./255., 226./255., 97./255., 16./255., 0./255.};
2520  Double_t blue[9] = { 99./255., 142./255., 198./255., 201./255., 90./255., 22./255., 13./255., 8./255., 2./255.};
2521  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2522  }
2523  break;
2524 
2525  // Inverted Dark Body Radiator
2526  case 56:
2527  {
2528  Double_t red[9] = { 242./255., 234./255., 237./255., 230./255., 212./255., 156./255., 99./255., 45./255., 0./255.};
2529  Double_t green[9] = { 243./255., 238./255., 238./255., 168./255., 101./255., 45./255., 0./255., 0./255., 0./255.};
2530  Double_t blue[9] = { 230./255., 95./255., 11./255., 8./255., 9./255., 3./255., 1./255., 1./255., 0./255.};
2531  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2532  }
2533  break;
2534 
2535  // Bird
2536  case 57:
2537  {
2538  Double_t red[9] = { 0.2082, 0.0592, 0.0780, 0.0232, 0.1802, 0.5301, 0.8186, 0.9956, 0.9764};
2539  Double_t green[9] = { 0.1664, 0.3599, 0.5041, 0.6419, 0.7178, 0.7492, 0.7328, 0.7862, 0.9832};
2540  Double_t blue[9] = { 0.5293, 0.8684, 0.8385, 0.7914, 0.6425, 0.4662, 0.3499, 0.1968, 0.0539};
2541  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2542  }
2543  break;
2544 
2545  // Cubehelix
2546  case 58:
2547  {
2548  Double_t red[9] = { 0.0000, 0.0956, 0.0098, 0.2124, 0.6905, 0.9242, 0.7914, 0.7596, 1.0000};
2549  Double_t green[9] = { 0.0000, 0.1147, 0.3616, 0.5041, 0.4577, 0.4691, 0.6905, 0.9237, 1.0000};
2550  Double_t blue[9] = { 0.0000, 0.2669, 0.3121, 0.1318, 0.2236, 0.6741, 0.9882, 0.9593, 1.0000};
2551  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2552  }
2553  break;
2554 
2555  // Green Red Violet
2556  case 59:
2557  {
2558  Double_t red[9] = {13./255., 23./255., 25./255., 63./255., 76./255., 104./255., 137./255., 161./255., 206./255.};
2559  Double_t green[9] = {95./255., 67./255., 37./255., 21./255., 0./255., 12./255., 35./255., 52./255., 79./255.};
2560  Double_t blue[9] = { 4./255., 3./255., 2./255., 6./255., 11./255., 22./255., 49./255., 98./255., 208./255.};
2561  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2562  }
2563  break;
2564 
2565  // Blue Red Yellow
2566  case 60:
2567  {
2568  Double_t red[9] = {0./255., 61./255., 89./255., 122./255., 143./255., 160./255., 185./255., 204./255., 231./255.};
2569  Double_t green[9] = {0./255., 0./255., 0./255., 0./255., 14./255., 37./255., 72./255., 132./255., 235./255.};
2570  Double_t blue[9] = {0./255., 140./255., 224./255., 144./255., 4./255., 5./255., 6./255., 9./255., 13./255.};
2571  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2572  }
2573  break;
2574 
2575  // Ocean
2576  case 61:
2577  {
2578  Double_t red[9] = { 14./255., 7./255., 2./255., 0./255., 5./255., 11./255., 55./255., 131./255., 229./255.};
2579  Double_t green[9] = {105./255., 56./255., 26./255., 1./255., 42./255., 74./255., 131./255., 171./255., 229./255.};
2580  Double_t blue[9] = { 2./255., 21./255., 35./255., 60./255., 92./255., 113./255., 160./255., 185./255., 229./255.};
2581  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2582  }
2583  break;
2584 
2585  // Color Printable On Grey
2586  case 62:
2587  {
2588  Double_t red[9] = { 0./255., 0./255., 0./255., 70./255., 148./255., 231./255., 235./255., 237./255., 244./255.};
2589  Double_t green[9] = { 0./255., 0./255., 0./255., 0./255., 0./255., 69./255., 67./255., 216./255., 244./255.};
2590  Double_t blue[9] = { 0./255., 102./255., 228./255., 231./255., 177./255., 124./255., 137./255., 20./255., 244./255.};
2591  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2592  }
2593  break;
2594 
2595  // Alpine
2596  case 63:
2597  {
2598  Double_t red[9] = { 50./255., 56./255., 63./255., 68./255., 93./255., 121./255., 165./255., 192./255., 241./255.};
2599  Double_t green[9] = { 66./255., 81./255., 91./255., 96./255., 111./255., 128./255., 155./255., 189./255., 241./255.};
2600  Double_t blue[9] = { 97./255., 91./255., 75./255., 65./255., 77./255., 103./255., 143./255., 167./255., 217./255.};
2601  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2602  }
2603  break;
2604 
2605  // Aquamarine
2606  case 64:
2607  {
2608  Double_t red[9] = { 145./255., 166./255., 167./255., 156./255., 131./255., 114./255., 101./255., 112./255., 132./255.};
2609  Double_t green[9] = { 158./255., 178./255., 179./255., 181./255., 163./255., 154./255., 144./255., 152./255., 159./255.};
2610  Double_t blue[9] = { 190./255., 199./255., 201./255., 192./255., 176./255., 169./255., 160./255., 166./255., 190./255.};
2611  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2612  }
2613  break;
2614 
2615  // Army
2616  case 65:
2617  {
2618  Double_t red[9] = { 93./255., 91./255., 99./255., 108./255., 130./255., 125./255., 132./255., 155./255., 174./255.};
2619  Double_t green[9] = { 126./255., 124./255., 128./255., 129./255., 131./255., 121./255., 119./255., 153./255., 173./255.};
2620  Double_t blue[9] = { 103./255., 94./255., 87./255., 85./255., 80./255., 85./255., 107./255., 120./255., 146./255.};
2621  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2622  }
2623  break;
2624 
2625  // Atlantic
2626  case 66:
2627  {
2628  Double_t red[9] = { 24./255., 40./255., 69./255., 90./255., 104./255., 114./255., 120./255., 132./255., 103./255.};
2629  Double_t green[9] = { 29./255., 52./255., 94./255., 127./255., 150./255., 162./255., 159./255., 151./255., 101./255.};
2630  Double_t blue[9] = { 29./255., 52./255., 96./255., 132./255., 162./255., 181./255., 184./255., 186./255., 131./255.};
2631  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2632  }
2633  break;
2634 
2635  // Aurora
2636  case 67:
2637  {
2638  Double_t red[9] = { 46./255., 38./255., 61./255., 92./255., 113./255., 121./255., 132./255., 150./255., 191./255.};
2639  Double_t green[9] = { 46./255., 36./255., 40./255., 69./255., 110./255., 135./255., 131./255., 92./255., 34./255.};
2640  Double_t blue[9] = { 46./255., 80./255., 74./255., 70./255., 81./255., 105./255., 165./255., 211./255., 225./255.};
2641  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2642  }
2643  break;
2644 
2645  // Avocado
2646  case 68:
2647  {
2648  Double_t red[9] = { 0./255., 4./255., 12./255., 30./255., 52./255., 101./255., 142./255., 190./255., 237./255.};
2649  Double_t green[9] = { 0./255., 40./255., 86./255., 121./255., 140./255., 172./255., 187./255., 213./255., 240./255.};
2650  Double_t blue[9] = { 0./255., 9./255., 14./255., 18./255., 21./255., 23./255., 27./255., 35./255., 101./255.};
2651  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2652  }
2653  break;
2654 
2655  // Beach
2656  case 69:
2657  {
2658  Double_t red[9] = { 198./255., 206./255., 206./255., 211./255., 198./255., 181./255., 161./255., 171./255., 244./255.};
2659  Double_t green[9] = { 103./255., 133./255., 150./255., 172./255., 178./255., 174./255., 163./255., 175./255., 244./255.};
2660  Double_t blue[9] = { 49./255., 54./255., 55./255., 66./255., 91./255., 130./255., 184./255., 224./255., 244./255.};
2661  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2662  }
2663  break;
2664 
2665  // Black Body
2666  case 70:
2667  {
2668  Double_t red[9] = { 243./255., 243./255., 240./255., 240./255., 241./255., 239./255., 186./255., 151./255., 129./255.};
2669  Double_t green[9] = { 0./255., 46./255., 99./255., 149./255., 194./255., 220./255., 183./255., 166./255., 147./255.};
2670  Double_t blue[9] = { 6./255., 8./255., 36./255., 91./255., 169./255., 235./255., 246./255., 240./255., 233./255.};
2671  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2672  }
2673  break;
2674 
2675  // Blue Green Yellow
2676  case 71:
2677  {
2678  Double_t red[9] = { 22./255., 19./255., 19./255., 25./255., 35./255., 53./255., 88./255., 139./255., 210./255.};
2679  Double_t green[9] = { 0./255., 32./255., 69./255., 108./255., 135./255., 159./255., 183./255., 198./255., 215./255.};
2680  Double_t blue[9] = { 77./255., 96./255., 110./255., 116./255., 110./255., 100./255., 90./255., 78./255., 70./255.};
2681  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2682  }
2683  break;
2684 
2685  // Brown Cyan
2686  case 72:
2687  {
2688  Double_t red[9] = { 68./255., 116./255., 165./255., 182./255., 189./255., 180./255., 145./255., 111./255., 71./255.};
2689  Double_t green[9] = { 37./255., 82./255., 135./255., 178./255., 204./255., 225./255., 221./255., 202./255., 147./255.};
2690  Double_t blue[9] = { 16./255., 55./255., 105./255., 147./255., 196./255., 226./255., 232./255., 224./255., 178./255.};
2691  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2692  }
2693  break;
2694 
2695  // CMYK
2696  case 73:
2697  {
2698  Double_t red[9] = { 61./255., 99./255., 136./255., 181./255., 213./255., 225./255., 198./255., 136./255., 24./255.};
2699  Double_t green[9] = { 149./255., 140./255., 96./255., 83./255., 132./255., 178./255., 190./255., 135./255., 22./255.};
2700  Double_t blue[9] = { 214./255., 203./255., 168./255., 135./255., 110./255., 100./255., 111./255., 113./255., 22./255.};
2701  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2702  }
2703  break;
2704 
2705  // Candy
2706  case 74:
2707  {
2708  Double_t red[9] = { 76./255., 120./255., 156./255., 183./255., 197./255., 180./255., 162./255., 154./255., 140./255.};
2709  Double_t green[9] = { 34./255., 35./255., 42./255., 69./255., 102./255., 137./255., 164./255., 188./255., 197./255.};
2710  Double_t blue[9] = { 64./255., 69./255., 78./255., 105./255., 142./255., 177./255., 205./255., 217./255., 198./255.};
2711  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2712  }
2713  break;
2714 
2715  // Cherry
2716  case 75:
2717  {
2718  Double_t red[9] = { 37./255., 102./255., 157./255., 188./255., 196./255., 214./255., 223./255., 235./255., 251./255.};
2719  Double_t green[9] = { 37./255., 29./255., 25./255., 37./255., 67./255., 91./255., 132./255., 185./255., 251./255.};
2720  Double_t blue[9] = { 37./255., 32./255., 33./255., 45./255., 66./255., 98./255., 137./255., 187./255., 251./255.};
2721  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2722  }
2723  break;
2724 
2725  // Coffee
2726  case 76:
2727  {
2728  Double_t red[9] = { 79./255., 100./255., 119./255., 137./255., 153./255., 172./255., 192./255., 205./255., 250./255.};
2729  Double_t green[9] = { 63./255., 79./255., 93./255., 103./255., 115./255., 135./255., 167./255., 196./255., 250./255.};
2730  Double_t blue[9] = { 51./255., 59./255., 66./255., 61./255., 62./255., 70./255., 110./255., 160./255., 250./255.};
2731  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2732  }
2733  break;
2734 
2735  // Dark Rain Bow
2736  case 77:
2737  {
2738  Double_t red[9] = { 43./255., 44./255., 50./255., 66./255., 125./255., 172./255., 178./255., 155./255., 157./255.};
2739  Double_t green[9] = { 63./255., 63./255., 85./255., 101./255., 138./255., 163./255., 122./255., 51./255., 39./255.};
2740  Double_t blue[9] = { 121./255., 101./255., 58./255., 44./255., 47./255., 55./255., 57./255., 44./255., 43./255.};
2741  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2742  }
2743  break;
2744 
2745  // Dark Terrain
2746  case 78:
2747  {
2748  Double_t red[9] = { 0./255., 41./255., 62./255., 79./255., 90./255., 87./255., 99./255., 140./255., 228./255.};
2749  Double_t green[9] = { 0./255., 57./255., 81./255., 93./255., 85./255., 70./255., 71./255., 125./255., 228./255.};
2750  Double_t blue[9] = { 95./255., 91./255., 91./255., 82./255., 60./255., 43./255., 44./255., 112./255., 228./255.};
2751  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2752  }
2753  break;
2754 
2755  // Fall
2756  case 79:
2757  {
2758  Double_t red[9] = { 49./255., 59./255., 72./255., 88./255., 114./255., 141./255., 176./255., 205./255., 222./255.};
2759  Double_t green[9] = { 78./255., 72./255., 66./255., 57./255., 59./255., 75./255., 106./255., 142./255., 173./255.};
2760  Double_t blue[9] = { 78./255., 55./255., 46./255., 40./255., 39./255., 39./255., 40./255., 41./255., 47./255.};
2761  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2762  }
2763  break;
2764 
2765  // Fruit Punch
2766  case 80:
2767  {
2768  Double_t red[9] = { 243./255., 222./255., 201./255., 185./255., 165./255., 158./255., 166./255., 187./255., 219./255.};
2769  Double_t green[9] = { 94./255., 108./255., 132./255., 135./255., 125./255., 96./255., 68./255., 51./255., 61./255.};
2770  Double_t blue[9] = { 7./255., 9./255., 12./255., 19./255., 45./255., 89./255., 118./255., 146./255., 118./255.};
2771  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2772  }
2773  break;
2774 
2775  // Fuchsia
2776  case 81:
2777  {
2778  Double_t red[9] = { 19./255., 44./255., 74./255., 105./255., 137./255., 166./255., 194./255., 206./255., 220./255.};
2779  Double_t green[9] = { 19./255., 28./255., 40./255., 55./255., 82./255., 110./255., 159./255., 181./255., 220./255.};
2780  Double_t blue[9] = { 19./255., 42./255., 68./255., 96./255., 129./255., 157./255., 188./255., 203./255., 220./255.};
2781  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2782  }
2783  break;
2784 
2785  // Grey Yellow
2786  case 82:
2787  {
2788  Double_t red[9] = { 33./255., 44./255., 70./255., 99./255., 140./255., 165./255., 199./255., 211./255., 216./255.};
2789  Double_t green[9] = { 38./255., 50./255., 76./255., 105./255., 140./255., 165./255., 191./255., 189./255., 167./255.};
2790  Double_t blue[9] = { 55./255., 67./255., 97./255., 124./255., 140./255., 166./255., 163./255., 129./255., 52./255.};
2791  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2792  }
2793  break;
2794 
2795  // Green Brown Terrain
2796  case 83:
2797  {
2798  Double_t red[9] = { 0./255., 33./255., 73./255., 124./255., 136./255., 152./255., 159./255., 171./255., 223./255.};
2799  Double_t green[9] = { 0./255., 43./255., 92./255., 124./255., 134./255., 126./255., 121./255., 144./255., 223./255.};
2800  Double_t blue[9] = { 0./255., 43./255., 68./255., 76./255., 73./255., 64./255., 72./255., 114./255., 223./255.};
2801  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2802  }
2803  break;
2804 
2805  // Green Pink
2806  case 84:
2807  {
2808  Double_t red[9] = { 5./255., 18./255., 45./255., 124./255., 193./255., 223./255., 205./255., 128./255., 49./255.};
2809  Double_t green[9] = { 48./255., 134./255., 207./255., 230./255., 193./255., 113./255., 28./255., 0./255., 7./255.};
2810  Double_t blue[9] = { 6./255., 15./255., 41./255., 121./255., 193./255., 226./255., 208./255., 130./255., 49./255.};
2811  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2812  }
2813  break;
2814 
2815  // Island
2816  case 85:
2817  {
2818  Double_t red[9] = { 180./255., 106./255., 104./255., 135./255., 164./255., 188./255., 189./255., 165./255., 144./255.};
2819  Double_t green[9] = { 72./255., 126./255., 154./255., 184./255., 198./255., 207./255., 205./255., 190./255., 179./255.};
2820  Double_t blue[9] = { 41./255., 120./255., 158./255., 188./255., 194./255., 181./255., 145./255., 100./255., 62./255.};
2821  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2822  }
2823  break;
2824 
2825  // Lake
2826  case 86:
2827  {
2828  Double_t red[9] = { 57./255., 72./255., 94./255., 117./255., 136./255., 154./255., 174./255., 192./255., 215./255.};
2829  Double_t green[9] = { 0./255., 33./255., 68./255., 109./255., 140./255., 171./255., 192./255., 196./255., 209./255.};
2830  Double_t blue[9] = { 116./255., 137./255., 173./255., 201./255., 200./255., 201./255., 203./255., 190./255., 187./255.};
2831  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2832  }
2833  break;
2834 
2835  // Light Temperature
2836  case 87:
2837  {
2838  Double_t red[9] = { 31./255., 71./255., 123./255., 160./255., 210./255., 222./255., 214./255., 199./255., 183./255.};
2839  Double_t green[9] = { 40./255., 117./255., 171./255., 211./255., 231./255., 220./255., 190./255., 132./255., 65./255.};
2840  Double_t blue[9] = { 234./255., 214./255., 228./255., 222./255., 210./255., 160./255., 105./255., 60./255., 34./255.};
2841  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2842  }
2843  break;
2844 
2845  // Light Terrain
2846  case 88:
2847  {
2848  Double_t red[9] = { 123./255., 108./255., 109./255., 126./255., 154./255., 172./255., 188./255., 196./255., 218./255.};
2849  Double_t green[9] = { 184./255., 138./255., 130./255., 133./255., 154./255., 175./255., 188./255., 196./255., 218./255.};
2850  Double_t blue[9] = { 208./255., 130./255., 109./255., 99./255., 110./255., 122./255., 150./255., 171./255., 218./255.};
2851  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2852  }
2853  break;
2854 
2855  // Mint
2856  case 89:
2857  {
2858  Double_t red[9] = { 105./255., 106./255., 122./255., 143./255., 159./255., 172./255., 176./255., 181./255., 207./255.};
2859  Double_t green[9] = { 252./255., 197./255., 194./255., 187./255., 174./255., 162./255., 153./255., 136./255., 125./255.};
2860  Double_t blue[9] = { 146./255., 133./255., 144./255., 155./255., 163./255., 167./255., 166./255., 162./255., 174./255.};
2861  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2862  }
2863  break;
2864 
2865  // Neon
2866  case 90:
2867  {
2868  Double_t red[9] = { 171./255., 141./255., 145./255., 152./255., 154./255., 159./255., 163./255., 158./255., 177./255.};
2869  Double_t green[9] = { 236./255., 143./255., 100./255., 63./255., 53./255., 55./255., 44./255., 31./255., 6./255.};
2870  Double_t blue[9] = { 59./255., 48./255., 46./255., 44./255., 42./255., 54./255., 82./255., 112./255., 179./255.};
2871  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2872  }
2873  break;
2874 
2875  // Pastel
2876  case 91:
2877  {
2878  Double_t red[9] = { 180./255., 190./255., 209./255., 223./255., 204./255., 228./255., 205./255., 152./255., 91./255.};
2879  Double_t green[9] = { 93./255., 125./255., 147./255., 172./255., 181./255., 224./255., 233./255., 198./255., 158./255.};
2880  Double_t blue[9] = { 236./255., 218./255., 160./255., 133./255., 114./255., 132./255., 162./255., 220./255., 218./255.};
2881  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2882  }
2883  break;
2884 
2885  // Pearl
2886  case 92:
2887  {
2888  Double_t red[9] = { 225./255., 183./255., 162./255., 135./255., 115./255., 111./255., 119./255., 145./255., 211./255.};
2889  Double_t green[9] = { 205./255., 177./255., 166./255., 135./255., 124./255., 117./255., 117./255., 132./255., 172./255.};
2890  Double_t blue[9] = { 186./255., 165./255., 155./255., 135./255., 126./255., 130./255., 150./255., 178./255., 226./255.};
2891  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2892  }
2893  break;
2894 
2895  // Pigeon
2896  case 93:
2897  {
2898  Double_t red[9] = { 39./255., 43./255., 59./255., 63./255., 80./255., 116./255., 153./255., 177./255., 223./255.};
2899  Double_t green[9] = { 39./255., 43./255., 59./255., 74./255., 91./255., 114./255., 139./255., 165./255., 223./255.};
2900  Double_t blue[9] = { 39./255., 50./255., 59./255., 70./255., 85./255., 115./255., 151./255., 176./255., 223./255.};
2901  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2902  }
2903  break;
2904 
2905  // Plum
2906  case 94:
2907  {
2908  Double_t red[9] = { 0./255., 38./255., 60./255., 76./255., 84./255., 89./255., 101./255., 128./255., 204./255.};
2909  Double_t green[9] = { 0./255., 10./255., 15./255., 23./255., 35./255., 57./255., 83./255., 123./255., 199./255.};
2910  Double_t blue[9] = { 0./255., 11./255., 22./255., 40./255., 63./255., 86./255., 97./255., 94./255., 85./255.};
2911  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2912  }
2913  break;
2914 
2915  // Red Blue
2916  case 95:
2917  {
2918  Double_t red[9] = { 94./255., 112./255., 141./255., 165./255., 167./255., 140./255., 91./255., 49./255., 27./255.};
2919  Double_t green[9] = { 27./255., 46./255., 88./255., 135./255., 166./255., 161./255., 135./255., 97./255., 58./255.};
2920  Double_t blue[9] = { 42./255., 52./255., 81./255., 106./255., 139./255., 158./255., 155./255., 137./255., 116./255.};
2921  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2922  }
2923  break;
2924 
2925  // Rose
2926  case 96:
2927  {
2928  Double_t red[9] = { 30./255., 49./255., 79./255., 117./255., 135./255., 151./255., 146./255., 138./255., 147./255.};
2929  Double_t green[9] = { 63./255., 60./255., 72./255., 90./255., 94./255., 94./255., 68./255., 46./255., 16./255.};
2930  Double_t blue[9] = { 18./255., 28./255., 41./255., 56./255., 62./255., 63./255., 50./255., 36./255., 21./255.};
2931  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2932  }
2933  break;
2934 
2935  // Rust
2936  case 97:
2937  {
2938  Double_t red[9] = { 0./255., 30./255., 63./255., 101./255., 143./255., 152./255., 169./255., 187./255., 230./255.};
2939  Double_t green[9] = { 0./255., 14./255., 28./255., 42./255., 58./255., 61./255., 67./255., 74./255., 91./255.};
2940  Double_t blue[9] = { 39./255., 26./255., 21./255., 18./255., 15./255., 14./255., 14./255., 13./255., 13./255.};
2941  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2942  }
2943  break;
2944 
2945  // Sandy Terrain
2946  case 98:
2947  {
2948  Double_t red[9] = { 149./255., 140./255., 164./255., 179./255., 182./255., 181./255., 131./255., 87./255., 61./255.};
2949  Double_t green[9] = { 62./255., 70./255., 107./255., 136./255., 144./255., 138./255., 117./255., 87./255., 74./255.};
2950  Double_t blue[9] = { 40./255., 38./255., 45./255., 49./255., 49./255., 49./255., 38./255., 32./255., 34./255.};
2951  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2952  }
2953  break;
2954 
2955  // Sienna
2956  case 99:
2957  {
2958  Double_t red[9] = { 99./255., 112./255., 148./255., 165./255., 179./255., 182./255., 183./255., 183./255., 208./255.};
2959  Double_t green[9] = { 39./255., 40./255., 57./255., 79./255., 104./255., 127./255., 148./255., 161./255., 198./255.};
2960  Double_t blue[9] = { 15./255., 16./255., 18./255., 33./255., 51./255., 79./255., 103./255., 129./255., 177./255.};
2961  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2962  }
2963  break;
2964 
2965  // Solar
2966  case 100:
2967  {
2968  Double_t red[9] = { 99./255., 116./255., 154./255., 174./255., 200./255., 196./255., 201./255., 201./255., 230./255.};
2969  Double_t green[9] = { 0./255., 0./255., 8./255., 32./255., 58./255., 83./255., 119./255., 136./255., 173./255.};
2970  Double_t blue[9] = { 5./255., 6./255., 7./255., 9./255., 9./255., 14./255., 17./255., 19./255., 24./255.};
2971  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2972  }
2973  break;
2974 
2975  // South West
2976  case 101:
2977  {
2978  Double_t red[9] = { 82./255., 106./255., 126./255., 141./255., 155./255., 163./255., 142./255., 107./255., 66./255.};
2979  Double_t green[9] = { 62./255., 44./255., 69./255., 107./255., 135./255., 152./255., 149./255., 132./255., 119./255.};
2980  Double_t blue[9] = { 39./255., 25./255., 31./255., 60./255., 73./255., 68./255., 49./255., 72./255., 188./255.};
2981  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2982  }
2983  break;
2984 
2985  // Starry Night
2986  case 102:
2987  {
2988  Double_t red[9] = { 18./255., 29./255., 44./255., 72./255., 116./255., 158./255., 184./255., 208./255., 221./255.};
2989  Double_t green[9] = { 27./255., 46./255., 71./255., 105./255., 146./255., 177./255., 189./255., 190./255., 183./255.};
2990  Double_t blue[9] = { 39./255., 55./255., 80./255., 108./255., 130./255., 133./255., 124./255., 100./255., 76./255.};
2991  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
2992  }
2993  break;
2994 
2995  // Sunset
2996  case 103:
2997  {
2998  Double_t red[9] = { 0./255., 48./255., 119./255., 173./255., 212./255., 224./255., 228./255., 228./255., 245./255.};
2999  Double_t green[9] = { 0./255., 13./255., 30./255., 47./255., 79./255., 127./255., 167./255., 205./255., 245./255.};
3000  Double_t blue[9] = { 0./255., 68./255., 75./255., 43./255., 16./255., 22./255., 55./255., 128./255., 245./255.};
3001  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
3002  }
3003  break;
3004 
3005  // Temperature Map
3006  case 104:
3007  {
3008  Double_t red[9] = { 34./255., 70./255., 129./255., 187./255., 225./255., 226./255., 216./255., 193./255., 179./255.};
3009  Double_t green[9] = { 48./255., 91./255., 147./255., 194./255., 226./255., 229./255., 196./255., 110./255., 12./255.};
3010  Double_t blue[9] = { 234./255., 212./255., 216./255., 224./255., 206./255., 110./255., 53./255., 40./255., 29./255.};
3011  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
3012  }
3013  break;
3014 
3015  // Thermometer
3016  case 105:
3017  {
3018  Double_t red[9] = { 30./255., 55./255., 103./255., 147./255., 174./255., 203./255., 188./255., 151./255., 105./255.};
3019  Double_t green[9] = { 0./255., 65./255., 138./255., 182./255., 187./255., 175./255., 121./255., 53./255., 9./255.};
3020  Double_t blue[9] = { 191./255., 202./255., 212./255., 208./255., 171./255., 140./255., 97./255., 57./255., 30./255.};
3021  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
3022  }
3023  break;
3024 
3025  // Valentine
3026  case 106:
3027  {
3028  Double_t red[9] = { 112./255., 97./255., 113./255., 125./255., 138./255., 159./255., 178./255., 188./255., 225./255.};
3029  Double_t green[9] = { 16./255., 17./255., 24./255., 37./255., 56./255., 81./255., 110./255., 136./255., 189./255.};
3030  Double_t blue[9] = { 38./255., 35./255., 46./255., 59./255., 78./255., 103./255., 130./255., 152./255., 201./255.};
3031  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
3032  }
3033  break;
3034 
3035  // Visible Spectrum
3036  case 107:
3037  {
3038  Double_t red[9] = { 18./255., 72./255., 5./255., 23./255., 29./255., 201./255., 200./255., 98./255., 29./255.};
3039  Double_t green[9] = { 0./255., 0./255., 43./255., 167./255., 211./255., 117./255., 0./255., 0./255., 0./255.};
3040  Double_t blue[9] = { 51./255., 203./255., 177./255., 26./255., 10./255., 9./255., 8./255., 3./255., 0./255.};
3041  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
3042  }
3043  break;
3044 
3045  // Water Melon
3046  case 108:
3047  {
3048  Double_t red[9] = { 19./255., 42./255., 64./255., 88./255., 118./255., 147./255., 175./255., 187./255., 205./255.};
3049  Double_t green[9] = { 19./255., 55./255., 89./255., 125./255., 154./255., 169./255., 161./255., 129./255., 70./255.};
3050  Double_t blue[9] = { 19./255., 32./255., 47./255., 70./255., 100./255., 128./255., 145./255., 130./255., 75./255.};
3051  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
3052  }
3053  break;
3054 
3055  // Cool
3056  case 109:
3057  {
3058  Double_t red[9] = { 33./255., 31./255., 42./255., 68./255., 86./255., 111./255., 141./255., 172./255., 227./255.};
3059  Double_t green[9] = { 255./255., 175./255., 145./255., 106./255., 88./255., 55./255., 15./255., 0./255., 0./255.};
3060  Double_t blue[9] = { 255./255., 205./255., 202./255., 203./255., 208./255., 205./255., 203./255., 206./255., 231./255.};
3061  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
3062  }
3063  break;
3064 
3065  // Copper
3066  case 110:
3067  {
3068  Double_t red[9] = { 0./255., 25./255., 50./255., 79./255., 110./255., 145./255., 181./255., 201./255., 254./255.};
3069  Double_t green[9] = { 0./255., 16./255., 30./255., 46./255., 63./255., 82./255., 101./255., 124./255., 179./255.};
3070  Double_t blue[9] = { 0./255., 12./255., 21./255., 29./255., 39./255., 49./255., 61./255., 74./255., 103./255.};
3071  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
3072  }
3073  break;
3074 
3075  // Gist Earth
3076  case 111:
3077  {
3078  Double_t red[9] = { 0./255., 13./255., 30./255., 44./255., 72./255., 120./255., 156./255., 200./255., 247./255.};
3079  Double_t green[9] = { 0./255., 36./255., 84./255., 117./255., 141./255., 153./255., 151./255., 158./255., 247./255.};
3080  Double_t blue[9] = { 0./255., 94./255., 100./255., 82./255., 56./255., 66./255., 76./255., 131./255., 247./255.};
3081  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
3082  }
3083  break;
3084 
3085  // Viridis
3086  case 112:
3087  {
3088  Double_t red[9] = { 26./255., 51./255., 43./255., 33./255., 28./255., 35./255., 74./255., 144./255., 246./255.};
3089  Double_t green[9] = { 9./255., 24./255., 55./255., 87./255., 118./255., 150./255., 180./255., 200./255., 222./255.};
3090  Double_t blue[9] = { 30./255., 96./255., 112./255., 114./255., 112./255., 101./255., 72./255., 35./255., 0./255.};
3091  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
3092  }
3093  break;
3094 
3095  // Cividis
3096  case 113:
3097  {
3098  Double_t red[9] = { 0./255., 5./255., 65./255., 97./255., 124./255., 156./255., 189./255., 224./255., 255./255.};
3099  Double_t green[9] = { 32./255., 54./255., 77./255., 100./255., 123./255., 148./255., 175./255., 203./255., 234./255.};
3100  Double_t blue[9] = { 77./255., 110./255., 107./255., 111./255., 120./255., 119./255., 111./255., 94./255., 70./255.};
3101  Idx = TColor::CreateGradientColorTable(9, stops, red, green, blue, 255, alpha);
3102  }
3103  break;
3104 
3105  default:
3106  ::Error("SetPalette", "Unknown palette number %d", ncolors);
3107  return;
3108  }
3109  paletteType = ncolors;
3110  if (Idx>0) fgPalettesList.fArray[paletteType-51] = (Double_t)Idx;
3111  else fgPalettesList.fArray[paletteType-51] = 0.;
3112  if (alpha > 0.) fgPalettesList.fArray[paletteType-51] += alpha/10.;
3113  return;
3114  }
3115 
3116  // set user defined palette
3117  if (colors) {
3118  fgPalette.Set(ncolors);
3119  for (i=0;i<ncolors;i++) fgPalette.fArray[i] = colors[i];
3120  } else {
3121  fgPalette.Set(TMath::Min(50,ncolors));
3122  for (i=0;i<TMath::Min(50,ncolors);i++) fgPalette.fArray[i] = palette[i];
3123  }
3124  paletteType = 3;
3125 }
3126 
3127 
3128 ////////////////////////////////////////////////////////////////////////////////
3129 /// Invert the current color palette.
3130 /// The top of the palette becomes the bottom and vice versa.
3131 
3133 {
3134  std::reverse(fgPalette.fArray, fgPalette.fArray + fgPalette.GetSize());
3135 }
UShort_t fBlue
Definition: GuiTypes.h:313
static void HSV2RGB(Float_t h, Float_t s, Float_t v, Float_t &r, Float_t &g, Float_t &b)
Static method to compute RGB from HSV:
Definition: TColor.cxx:1516
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
static Int_t gLastDefinedColors
Previous number of defined colors.
Definition: TColor.cxx:45
#define fgPalettesList
Definition: TColor.cxx:49
Float_t fBlue
Fraction of Blue.
Definition: TColor.h:25
virtual void SetAlpha(Float_t a)
Definition: TColor.h:67
An array of TObjects.
Definition: TObjArray.h:37
unsigned long ULong_t
Definition: CPyCppyy.h:51
static void SetColorThreshold(Float_t t)
This method specifies the color threshold used by GetColor to retrieve a color.
Definition: TColor.cxx:1837
static Bool_t IsGrayscale()
Return whether all colors return grayscale values.
Definition: TColor.cxx:2155
Float_t GetRed() const
Definition: TColor.h:57
#define fgPalette
Definition: TColor.cxx:48
static constexpr double gray
Float_t fHue
Hue.
Definition: TColor.h:26
#define fgGrayscaleMode
Definition: TColor.cxx:47
float Float_t
Definition: RtypesCore.h:55
unsigned int UInt_t
Definition: CPyCppyy.h:44
const char Option_t
Definition: RtypesCore.h:64
#define s0(x)
Definition: RSha256.hxx:90
Definition: Rtypes.h:64
static const TArrayI & GetPalette()
Static function returning the current active palette.
Definition: TColor.cxx:1401
#define g(i)
Definition: RSha256.hxx:105
static Int_t GetColorDark(Int_t color)
Static function: Returns the dark color number corresponding to n If the TColor object does not exist...
Definition: TColor.cxx:1944
Float_t GetAlpha() const
Definition: TColor.h:63
unsigned short UShort_t
Definition: RtypesCore.h:38
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition: TNamed.cxx:140
Definition: Rtypes.h:63
static void SaveColor(std::ostream &out, Int_t ci)
Save a color with index > 228 as a C++ statement(s) on output stream out.
Definition: TColor.cxx:2117
int Int_t
Definition: CPyCppyy.h:43
TColor & operator=(const TColor &color)
Definition: TColor.cxx:1071
static Int_t GetFreeColorIndex()
Static function: Returns a free color index which can be used to define a user custom color...
Definition: TColor.cxx:2001
Definition: Rtypes.h:63
const Mask_t kDoRed
Definition: GuiTypes.h:318
Definition: Rtypes.h:65
static void RGB2HSV(Float_t r, Float_t g, Float_t b, Float_t &h, Float_t &s, Float_t &v)
Static method to compute HSV from RGB.
Definition: TColor.cxx:1648
const Mask_t kDoGreen
Definition: GuiTypes.h:319
static Int_t GetNumberOfColors()
Static function returning number of colors in the color palette.
Definition: TColor.cxx:1409
#define gROOT
Definition: TROOT.h:405
Basic string class.
Definition: TString.h:131
Float_t GetLight() const
Definition: TColor.h:61
#define f(i)
Definition: RSha256.hxx:104
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:180
bool Bool_t
Definition: RtypesCore.h:61
static void HLStoRGB(Float_t h, Float_t l, Float_t s, Float_t &r, Float_t &g, Float_t &b)
Definition: TColor.h:73
Definition: Rtypes.h:64
Definition: Rtypes.h:64
R__EXTERN TApplication * gApplication
Definition: TApplication.h:166
static Bool_t DefinedColors()
Static function returning kTRUE if some new colors have been defined after initialisation or since th...
Definition: TColor.cxx:1419
Float_t GetBlue() const
Definition: TColor.h:59
Short_t Abs(Short_t d)
Definition: TMathBase.h:120
Array of integers (32 bits per element).
Definition: TArrayI.h:27
Float_t fAlpha
Alpha (transparency)
Definition: TColor.h:29
static void SetPalette(Int_t ncolors, Int_t *colors, Float_t alpha=1.)
Static function.
Definition: TColor.cxx:2413
static const char * PixelAsHexString(ULong_t pixel)
Convert machine dependent pixel value (obtained via RGB2Pixel or via Number2Pixel() or via TColor::Ge...
Definition: TColor.cxx:2105
UShort_t fRed
Definition: GuiTypes.h:311
static Int_t GetColorPalette(Int_t i)
Static function returning the color number i in current palette.
Definition: TColor.cxx:1389
void Class()
Definition: Class.C:29
Definition: Rtypes.h:65
Float_t GetGreen() const
Definition: TColor.h:58
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:29
Float_t fLight
Light.
Definition: TColor.h:27
Float_t fGreen
Fraction of Green.
Definition: TColor.h:24
static constexpr double s
static ULong_t RGB2Pixel(Int_t r, Int_t g, Int_t b)
Convert r,g,b to graphics system dependent pixel value.
Definition: TColor.cxx:2050
Definition: Rtypes.h:65
virtual void Print(Option_t *option="") const
Dump this color with its attributes.
Definition: TColor.cxx:1581
virtual ~TColor()
Color destructor.
Definition: TColor.cxx:1057
static void RGB2HLS(Float_t r, Float_t g, Float_t b, Float_t &h, Float_t &l, Float_t &s)
Static method to compute HLS from RGB.
Definition: TColor.cxx:1590
void Copy(TObject &color) const
Copy this color to obj.
Definition: TColor.cxx:1229
static void RGBtoHLS(Float_t r, Float_t g, Float_t b, Float_t &h, Float_t &l, Float_t &s)
Definition: TColor.h:78
Definition: Rtypes.h:63
UShort_t fGreen
Definition: GuiTypes.h:312
static void HLS2RGB(Float_t h, Float_t l, Float_t s, Float_t &r, Float_t &g, Float_t &b)
Static method to compute RGB from HLS.
Definition: TColor.cxx:1450
static Int_t GetColorBright(Int_t color)
Static function: Returns the bright color number corresponding to n If the TColor object does not exi...
Definition: TColor.cxx:1912
static void InvertPalette()
Invert the current color palette.
Definition: TColor.cxx:3132
Float_t GetSaturation() const
Definition: TColor.h:62
virtual void AddAtAndExpand(TObject *obj, Int_t idx)
Add object at position idx.
Definition: TObjArray.cxx:235
static void NeedGraphicsLibs()
Static method.
ROOT::R::TRInterface & r
Definition: Object.C:4
Definition: Rtypes.h:65
Float_t fSaturation
Saturation.
Definition: TColor.h:28
static void CreateColorsRectangle(Int_t offset, const char *name, UChar_t *rgb)
Create the "rectangular" colors in the color wheel.
Definition: TColor.cxx:1281
auto * a
Definition: textangle.C:12
static Float_t gColorThreshold
Color threshold used by GetColor.
Definition: TColor.cxx:43
void Form(const char *fmt,...)
Formats a string using a printf style format descriptor.
Definition: TString.cxx:2289
ULong_t GetPixel() const
Return pixel value corresponding to this color.
Definition: TColor.cxx:1433
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:889
char * Form(const char *fmt,...)
double floor(double)
static Int_t GetColor(const char *hexcolor)
Static method returning color number for color specified by hex color string of form: "#rrggbb"...
Definition: TColor.cxx:1765
static Int_t gDefinedColors
Number of defined colors.
Definition: TColor.cxx:44
TColor()
Default constructor.
Definition: TColor.cxx:982
void Warning(const char *location, const char *msgfmt,...)
ULong_t fPixel
Definition: GuiTypes.h:310
#define gVirtualX
Definition: TVirtualX.h:338
#define h(i)
Definition: RSha256.hxx:106
static Int_t GetColorTransparent(Int_t color, Float_t a)
Static function: Returns the transparent color number corresponding to n.
Definition: TColor.cxx:1975
const Bool_t kFALSE
Definition: RtypesCore.h:90
Int_t fNumber
Color number identifier.
Definition: TColor.h:21
static void CreateColorWheel()
Static function steering the creation of all colors in the color wheel.
Definition: TColor.cxx:1301
void InitializeGraphics()
Initialize the graphics environment.
static ULong_t Number2Pixel(Int_t ci)
Static method that given a color index number, returns the corresponding pixel value.
Definition: TColor.cxx:2012
Color * colors
Definition: X3DBuffer.c:21
static Int_t gHighestColorIndex
Highest color index defined.
Definition: TColor.cxx:42
#define ClassImp(name)
Definition: Rtypes.h:361
double Double_t
Definition: RtypesCore.h:57
static void CreateColorsGray()
Create the Gray scale colors in the Color Wheel.
Definition: TColor.cxx:1245
virtual void SetRGB(Float_t r, Float_t g, Float_t b)
Initialize this color and its associated colors.
Definition: TColor.cxx:1701
The color creation and management class.
Definition: TColor.h:19
Definition: Rtypes.h:64
Array of doubles (64 bits per element).
Definition: TArrayD.h:27
Float_t fRed
Fraction of Red.
Definition: TColor.h:23
Mother of all ROOT objects.
Definition: TObject.h:37
const Mask_t kDoBlue
Definition: GuiTypes.h:320
UShort_t fMask
Definition: GuiTypes.h:314
virtual void Copy(TObject &named) const
Copy this to obj.
Definition: TNamed.cxx:94
auto * l
Definition: textangle.C:4
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:212
static void SetGrayscale(Bool_t set=kTRUE)
Set whether all colors should return grayscale values.
Definition: TColor.cxx:2163
static void CreateColorsCircle(Int_t offset, const char *name, UChar_t *rgb)
Create the "circle" colors in the color wheel.
Definition: TColor.cxx:1261
Definition: Rtypes.h:65
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
#define snprintf
Definition: civetweb.c:1540
Int_t GetNumber() const
Definition: TColor.h:55
virtual void ls(Option_t *option="") const
List this color with its attributes.
Definition: TColor.cxx:1572
#define c(i)
Definition: RSha256.hxx:101
static void InitializeColors()
Initialize colors used by the TCanvas based graphics (via TColor objects).
Definition: TColor.cxx:1082
Definition: Rtypes.h:64
unsigned char UChar_t
Definition: RtypesCore.h:36
void Allocate()
Make this color known to the graphics system.
Definition: TColor.cxx:1748
const char * AsHexString() const
Return color as hexadecimal string.
Definition: TColor.cxx:1208
float * q
Definition: THbookFile.cxx:87
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition: TNamed.cxx:164
const Bool_t kTRUE
Definition: RtypesCore.h:89
static Int_t CreateGradientColorTable(UInt_t Number, Double_t *Stops, Double_t *Red, Double_t *Green, Double_t *Blue, UInt_t NColors, Float_t alpha=1.)
Static function creating a color table with several connected linear gradients.
Definition: TColor.cxx:2224
const Int_t n
Definition: legend1.C:16
static Float_t HLStoRGB1(Float_t rn1, Float_t rn2, Float_t huei)
Static method. Auxiliary to HLS2RGB().
Definition: TColor.cxx:1475
char name[80]
Definition: TGX11.cxx:109
Definition: Rtypes.h:65
static void Pixel2RGB(ULong_t pixel, Int_t &r, Int_t &g, Int_t &b)
Convert machine dependent pixel value (obtained via RGB2Pixel or via Number2Pixel() or via TColor::Ge...
Definition: TColor.cxx:2088
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:875
Float_t GetHue() const
Definition: TColor.h:60
const char * Data() const
Definition: TString.h:364