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