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