source: trunk/MagicSoft/Mars/mhbase/MH.cc@ 8223

Last change on this file since 8223 was 8106, checked in by tbretz, 18 years ago
*** empty log message ***
File size: 44.3 KB
Line 
1/* ======================================================================== *\
2! $Name: not supported by cvs2svn $:$Id: MH.cc,v 1.31 2006-10-17 17:16:00 tbretz Exp $
3! --------------------------------------------------------------------------
4!
5! *
6! * This file is part of MARS, the MAGIC Analysis and Reconstruction
7! * Software. It is distributed to you in the hope that it can be a useful
8! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
9! * It is distributed WITHOUT ANY WARRANTY.
10! *
11! * Permission to use, copy, modify and distribute this software and its
12! * documentation for any purpose is hereby granted without fee,
13! * provided that the above copyright notice appear in all copies and
14! * that both that copyright notice and this permission notice appear
15! * in supporting documentation. It is provided "as is" without express
16! * or implied warranty.
17! *
18!
19!
20! Author(s): Thomas Bretz 07/2001 <mailto:tbretz@astro.uni-wuerzburg.de>
21!
22! Copyright: MAGIC Software Development, 2000-2002
23!
24!
25\* ======================================================================== */
26
27//////////////////////////////////////////////////////////////////////////////
28// //
29// MH //
30// //
31// This is a base tasks for mars histograms. It defines a common interface //
32// for filling the histograms with events (MH::Fill) which is used by a //
33// common 'filler' And a SetupFill member function which may be used //
34// by MFillH. The idea is: //
35// 1) If your Histogram can become filled by one single container //
36// (like MHHillas) you overload MH::Fill and it gets called with //
37// a pointer to the container with which it should be filled. //
38// //
39// 2) You histogram needs several containers to get filled. Than you //
40// have to overload MH::SetupFill and get the necessary objects from //
41// the parameter list. Use this objects in Fill to fill your //
42// histogram. //
43// //
44// If you want to create your own histogram class the new class must be //
45// derived from MH (instead of the base MParContainer) and you must //
46// the fill function of MH. This is the function which is called to fill //
47// the histogram(s) by the data of a corresponding parameter container. //
48// //
49// Remark: the static member function (eg MakeDefCanvas) can be called //
50// from everywhere using: MH::MakeDefCanvas(...) //
51// //
52//////////////////////////////////////////////////////////////////////////////
53
54#include "MH.h"
55
56#include <TH1.h>
57#include <TH2.h>
58#include <TH3.h>
59#include <TStyle.h> // TStyle::GetScreenFactor
60#include <TGaxis.h>
61#include <TCanvas.h>
62#include <TLegend.h>
63#include <TPaveStats.h>
64#include <TBaseClass.h>
65#if ROOT_VERSION_CODE > ROOT_VERSION(3,04,01)
66#include <THLimitsFinder.h>
67#endif
68
69#include "MLog.h"
70#include "MLogManip.h"
71
72#include "MParList.h"
73#include "MParContainer.h"
74
75#include "MBinning.h"
76
77#include "MArrayD.h"
78#include "MArrayF.h"
79
80ClassImp(MH);
81
82using namespace std;
83
84// --------------------------------------------------------------------------
85//
86// Default Constructor. It sets name and title only. Typically you won't
87// need to change this.
88//
89MH::MH(const char *name, const char *title)
90 : fSerialNumber(0), fNumExecutions(0)
91
92{
93 //
94 // set the name and title of this object
95 //
96 fName = name ? name : "MH";
97 fTitle = title ? title : "Base class for Mars histograms";
98}
99
100// --------------------------------------------------------------------------
101//
102// If you want to use the automatic filling of your derived class you
103// must overload this function. If it is not overloaded you cannot use
104// FillH with this class. The argument is a pointer to a container
105// in your paremeter list which is specified in the MFillH constructor.
106// If you are not going to use it you should at least add
107// Bool_t MH::Fill(const MParContainer *) { return kTRUE; }
108// to your class definition.
109//
110Bool_t MH::Fill(const MParContainer *par, const Stat_t w)
111{
112 *fLog << warn << GetDescriptor() << ": Fill not overloaded! Can't be used!" << endl;
113 return kFALSE;
114}
115
116// --------------------------------------------------------------------------
117//
118// This virtual function is ment as a generalized interface to retrieve
119// a pointer to a root histogram from the MH-derived class.
120//
121TH1 *MH::GetHistByName(const TString name) const
122{
123 *fLog << warn << GetDescriptor() << ": GetHistByName not overloaded! Can't be used!" << endl;
124 return NULL;
125}
126
127// --------------------------------------------------------------------------
128//
129// This is a function which should replace the creation of default
130// canvases like root does. Because this is inconvinient in some aspects.
131// need to change this.
132// You can specify a name for the default canvas and a title. Also
133// width and height can be given.
134// MakeDefCanvas looks for a canvas with the given name. If now name is
135// given the DefCanvasName of root is used. If no such canvas is existing
136// it is created and returned. If such a canvas already exists a new canvas
137// with a name plus anumber is created (the number is calculated by the
138// number of all existing canvases plus one)
139//
140// Normally the canvas size is scaled with gStyle->GetScreenFactor() so
141// that on all screens it looks like the same part of the screen.
142// To suppress this scaling use usescreenfactor=kFALSE. In this case
143// you specify directly the size of the embedded pad.
144//
145TCanvas *MH::MakeDefCanvas(TString name, const char *title,
146 UInt_t w, UInt_t h, Bool_t usescreenfactor)
147{
148 const TList *list = (TList*)gROOT->GetListOfCanvases();
149
150 if (name.IsNull())
151 name = gROOT->GetDefCanvasName();
152
153 if (list->FindObject(name))
154 name += Form(" <%d>", list->GetSize()+1);
155
156 if (!usescreenfactor)
157 {
158 const Float_t cx = gStyle->GetScreenFactor();
159 w += 4;
160 h += 28;
161 w = (int)(w/cx+1);
162 h = (int)(h/cx+1);
163 }
164
165 return new TCanvas(name, title, w, h);
166}
167
168// --------------------------------------------------------------------------
169//
170// This function works like MakeDefCanvas(name, title, w, h) but name
171// and title are retrieved from the given TObject.
172//
173// Normally the canvas size is scaled with gStyle->GetScreenFactor() so
174// that on all screens it looks like the same part of the screen.
175// To suppress this scaling use usescreenfactor=kFALSE. In this case
176// you specify directly the size of the embedded pad.
177//
178TCanvas *MH::MakeDefCanvas(const TObject *obj,
179 UInt_t w, UInt_t h, Bool_t usescreenfactor)
180{
181 if (!usescreenfactor)
182 {
183 const Float_t cx = gStyle->GetScreenFactor();
184 w += 4;
185 h += 28;
186 h = (int)(h/cx+1);
187 w = (int)(w/cx+1);
188 }
189
190 return MakeDefCanvas(obj->GetName(), obj->GetTitle(), w, h);
191}
192
193// --------------------------------------------------------------------------
194//
195// Search in gPad for all objects with the name name and remove all of them
196// (TList::Remove)
197//
198void MH::RemoveFromPad(const char *name)
199{
200 if (!gPad)
201 return;
202
203 TList *list = gPad->GetListOfPrimitives();
204 if (!list)
205 return;
206
207 TObject *obj = 0;
208 while ((obj = gPad->FindObject(name)))
209 list->Remove(obj);
210}
211
212// --------------------------------------------------------------------------
213//
214// Applies a given binning to a 1D-histogram
215//
216void MH::SetBinning(TH1 *h, const MBinning *binsx)
217{
218 //
219 // Another strange behaviour: TAxis::Set deletes the axis title!
220 //
221 TAxis &x = *h->GetXaxis();
222
223#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
224 TString xtitle = x.GetTitle();
225#endif
226
227 // All this is reset by TAxis::Set
228 const TAttAxis att(x);
229 const Bool_t tm(x.GetTimeDisplay());
230 const TString tf(x.GetTimeFormat());
231
232 //
233 // This is a necessary workaround if one wants to set
234 // non-equidistant bins after the initialization
235 // TH1D::fNcells must be set correctly.
236 //
237 h->SetBins(binsx->GetNumBins(), 0, 1);
238
239 //
240 // Set the binning of the current histogram to the binning
241 // in one of the two given histograms
242 //
243 x.Set(binsx->GetNumBins(), binsx->GetEdges());
244
245 // All this is reset by TAxis::Set
246 att.Copy(x);
247 x.SetTimeDisplay(tm);
248 x.SetTimeFormat(tf);
249
250#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
251 x.SetTitle(xtitle);
252#endif
253}
254
255// --------------------------------------------------------------------------
256//
257// Applies given binnings to the two axis of a 2D-histogram
258//
259void MH::SetBinning(TH2 *h, const MBinning *binsx, const MBinning *binsy)
260{
261 TAxis &x = *h->GetXaxis();
262 TAxis &y = *h->GetYaxis();
263
264 //
265 // Another strange behaviour: TAxis::Set deletes the axis title!
266 //
267#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
268 TString xtitle = x.GetTitle();
269 TString ytitle = y.GetTitle();
270#endif
271
272 // All this is reset by TAxis::Set
273 const TAttAxis attx(x);
274 const TAttAxis atty(y);
275 const Bool_t tmx(x.GetTimeDisplay());
276 const Bool_t tmy(y.GetTimeDisplay());
277 const TString tfx(x.GetTimeFormat());
278 const TString tfy(y.GetTimeFormat());
279
280 //
281 // This is a necessary workaround if one wants to set
282 // non-equidistant bins after the initialization
283 // TH1D::fNcells must be set correctly.
284 //
285 h->SetBins(binsx->GetNumBins(), 0, 1,
286 binsy->GetNumBins(), 0, 1);
287
288 //
289 // Set the binning of the current histogram to the binning
290 // in one of the two given histograms
291 //
292 x.Set(binsx->GetNumBins(), binsx->GetEdges());
293 y.Set(binsy->GetNumBins(), binsy->GetEdges());
294
295 // All this is reset by TAxis::Set
296 attx.Copy(x);
297 atty.Copy(y);
298 x.SetTimeDisplay(tmx);
299 y.SetTimeDisplay(tmy);
300 x.SetTimeFormat(tfx);
301 y.SetTimeFormat(tfy);
302
303#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
304 x.SetTitle(xtitle);
305 y.SetTitle(ytitle);
306#endif
307}
308
309// --------------------------------------------------------------------------
310//
311// Applies given binnings to the three axis of a 3D-histogram
312//
313void MH::SetBinning(TH3 *h, const MBinning *binsx, const MBinning *binsy, const MBinning *binsz)
314{
315 //
316 // Another strange behaviour: TAxis::Set deletes the axis title!
317 //
318 TAxis &x = *h->GetXaxis();
319 TAxis &y = *h->GetYaxis();
320 TAxis &z = *h->GetZaxis();
321
322#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
323 TString xtitle = x.GetTitle();
324 TString ytitle = y.GetTitle();
325 TString ztitle = z.GetTitle();
326#endif
327
328 // All this is reset by TAxis::Set
329 const TAttAxis attx(x);
330 const TAttAxis atty(y);
331 const TAttAxis attz(z);
332 const Bool_t tmx(x.GetTimeDisplay());
333 const Bool_t tmy(y.GetTimeDisplay());
334 const Bool_t tmz(z.GetTimeDisplay());
335 const TString tfx(x.GetTimeFormat());
336 const TString tfy(y.GetTimeFormat());
337 const TString tfz(z.GetTimeFormat());
338
339 //
340 // This is a necessary workaround if one wants to set
341 // non-equidistant bins after the initialization
342 // TH1D::fNcells must be set correctly.
343 //
344 h->SetBins(binsx->GetNumBins(), 0, 1,
345 binsy->GetNumBins(), 0, 1,
346 binsz->GetNumBins(), 0, 1);
347
348 //
349 // Set the binning of the current histogram to the binning
350 // in one of the two given histograms
351 //
352 x.Set(binsx->GetNumBins(), binsx->GetEdges());
353 y.Set(binsy->GetNumBins(), binsy->GetEdges());
354 z.Set(binsz->GetNumBins(), binsz->GetEdges());
355
356 // All this is reset by TAxis::Set
357 attx.Copy(x);
358 atty.Copy(y);
359 attz.Copy(z);
360 x.SetTimeDisplay(tmx);
361 y.SetTimeDisplay(tmy);
362 z.SetTimeDisplay(tmz);
363 x.SetTimeFormat(tfx);
364 y.SetTimeFormat(tfy);
365 z.SetTimeFormat(tfz);
366
367#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
368 x.SetTitle(xtitle);
369 y.SetTitle(ytitle);
370 z.SetTitle(ztitle);
371#endif
372}
373
374// --------------------------------------------------------------------------
375//
376// Applies given binning (the n+1 edges) to the axis of a 1D-histogram
377//
378void MH::SetBinning(TH1 *h, const TArrayD &binsx)
379{
380 MBinning bx;
381 bx.SetEdges(binsx);
382 SetBinning(h, &bx);
383}
384
385// --------------------------------------------------------------------------
386//
387// Applies given binning (the n+1 edges) to the two axis of a
388// 2D-histogram
389//
390void MH::SetBinning(TH2 *h, const TArrayD &binsx, const TArrayD &binsy)
391{
392 MBinning bx;
393 MBinning by;
394 bx.SetEdges(binsx);
395 by.SetEdges(binsy);
396 SetBinning(h, &bx, &by);
397}
398
399// --------------------------------------------------------------------------
400//
401// Applies given binning (the n+1 edges) to the three axis of a
402// 3D-histogram
403//
404void MH::SetBinning(TH3 *h, const TArrayD &binsx, const TArrayD &binsy, const TArrayD &binsz)
405{
406 MBinning bx;
407 MBinning by;
408 MBinning bz;
409 bx.SetEdges(binsx);
410 by.SetEdges(binsy);
411 bz.SetEdges(binsz);
412 SetBinning(h, &bx, &by, &bz);
413}
414
415// --------------------------------------------------------------------------
416//
417// Applies the binning of a TAxis (eg from a root histogram) to the axis
418// of a 1D-histogram
419//
420void MH::SetBinning(TH1 *h, const TAxis *binsx)
421{
422 const Int_t nx = binsx->GetNbins();
423
424 TArrayD bx(nx+1);
425 for (int i=0; i<nx; i++) bx[i] = binsx->GetBinLowEdge(i+1);
426 bx[nx] = binsx->GetXmax();
427
428 SetBinning(h, bx);
429}
430
431// --------------------------------------------------------------------------
432//
433// Applies the binnings of the TAxis' (eg from a root histogram) to the
434// two axis' of a 2D-histogram
435//
436void MH::SetBinning(TH2 *h, const TAxis *binsx, const TAxis *binsy)
437{
438 const Int_t nx = binsx->GetNbins();
439 const Int_t ny = binsy->GetNbins();
440
441 TArrayD bx(nx+1);
442 TArrayD by(ny+1);
443 for (int i=0; i<nx; i++) bx[i] = binsx->GetBinLowEdge(i+1);
444 for (int i=0; i<ny; i++) by[i] = binsy->GetBinLowEdge(i+1);
445 bx[nx] = binsx->GetXmax();
446 by[ny] = binsy->GetXmax();
447
448 SetBinning(h, bx, by);
449}
450
451// --------------------------------------------------------------------------
452//
453// Applies the binnings of the TAxis' (eg from a root histogram) to the
454// three axis' of a 3D-histogram
455//
456void MH::SetBinning(TH3 *h, const TAxis *binsx, const TAxis *binsy, const TAxis *binsz)
457{
458 const Int_t nx = binsx->GetNbins();
459 const Int_t ny = binsy->GetNbins();
460 const Int_t nz = binsz->GetNbins();
461
462 TArrayD bx(nx+1);
463 TArrayD by(ny+1);
464 TArrayD bz(nz+1);
465 for (int i=0; i<nx; i++) bx[i] = binsx->GetBinLowEdge(i+1);
466 for (int i=0; i<ny; i++) by[i] = binsy->GetBinLowEdge(i+1);
467 for (int i=0; i<nz; i++) bz[i] = binsz->GetBinLowEdge(i+1);
468 bx[nx] = binsx->GetXmax();
469 by[ny] = binsy->GetXmax();
470 bz[nz] = binsz->GetXmax();
471
472 SetBinning(h, bx, by, bz);
473}
474
475// --------------------------------------------------------------------------
476//
477// Applies the binnings of one root-histogram x to another one h
478// Both histograms must be of the same type: TH1, TH2 or TH3
479//
480void MH::SetBinning(TH1 *h, const TH1 *x)
481{
482 if (h->InheritsFrom(TH3::Class()) && x->InheritsFrom(TH3::Class()))
483 {
484 SetBinning((TH3*)h, ((TH1*)x)->GetXaxis(), ((TH1*)x)->GetYaxis(), ((TH1*)x)->GetZaxis());
485 return;
486 }
487 if (h->InheritsFrom(TH3::Class()) || x->InheritsFrom(TH3::Class()))
488 return;
489 if (h->InheritsFrom(TH2::Class()) && x->InheritsFrom(TH2::Class()))
490 {
491 SetBinning((TH2*)h, ((TH1*)x)->GetXaxis(), ((TH1*)x)->GetYaxis());
492 return;
493 }
494 if (h->InheritsFrom(TH2::Class()) || x->InheritsFrom(TH2::Class()))
495 return;
496 if (h->InheritsFrom(TH1::Class()) && x->InheritsFrom(TH1::Class()))
497 {
498 SetBinning(h, ((TH1*)x)->GetXaxis());
499 return;
500 }
501}
502
503void MH::RemoveFirstBin(TH1 &h)
504{
505 if (h.InheritsFrom(TH2::Class()) || h.InheritsFrom(TH3::Class()))
506 return;
507
508 const Bool_t haserr = h.GetSumw2N()>0;
509
510 const Int_t n0 = h.GetNbinsX();
511 if (n0<2)
512 return;
513
514 TArrayD val(n0-1);
515 TArrayD er(haserr ? n0-1 : 0);
516 for (int i=1; i<n0; i++)
517 {
518 val[i-1] = h.GetBinContent(i+1);
519 if (haserr)
520 er[i-1] = h.GetBinError(i+1);
521 }
522
523 MBinning bins;
524 bins.SetEdges(h, 'x');
525 bins.RemoveFirstEdge();
526 bins.Apply(h);
527
528 h.Reset();
529
530 for (int i=1; i<n0; i++)
531 {
532 h.SetBinContent(i, val[i-1]);
533 if (haserr)
534 h.SetBinError(i, er[i-1]);
535 }
536}
537
538// --------------------------------------------------------------------------
539//
540// Multiplies all entries in a TArrayD by a float f
541//
542void MH::ScaleArray(TArrayD &bins, Double_t f)
543{
544 for (int i=0; i<bins.GetSize(); i++)
545 bins[i] *= f;
546}
547
548// --------------------------------------------------------------------------
549//
550// Scales the binning of a TAxis by a float f
551//
552TArrayD MH::ScaleAxis(TAxis &axe, Double_t f)
553{
554 TArrayD arr(axe.GetNbins()+1);
555
556 for (int i=1; i<=axe.GetNbins()+1; i++)
557 arr[i-1] = axe.GetBinLowEdge(i);
558
559 ScaleArray(arr, f);
560
561 return arr;
562}
563
564// --------------------------------------------------------------------------
565//
566// Scales the binning of one, two or three axis of a histogram by a float f
567//
568void MH::ScaleAxis(TH1 *h, Double_t fx, Double_t fy, Double_t fz)
569{
570 if (h->InheritsFrom(TH3::Class()))
571 {
572 SetBinning((TH3*)h,
573 ScaleAxis(*h->GetXaxis(), fx),
574 ScaleAxis(*h->GetYaxis(), fy),
575 ScaleAxis(*h->GetZaxis(), fz));
576 return;
577 }
578
579 if (h->InheritsFrom(TH2::Class()))
580 {
581 SetBinning((TH2*)h,
582 ScaleAxis(*h->GetXaxis(), fx),
583 ScaleAxis(*h->GetYaxis(), fy));
584 return;
585 }
586
587 if (h->InheritsFrom(TH1::Class()))
588 SetBinning(h, ScaleAxis(*h->GetXaxis(), fx));
589}
590
591// --------------------------------------------------------------------------
592//
593// Tries to find a MBinning container with the name "Binning"+name
594// in the given parameter list. If it was found it is applied to the
595// given histogram. This is only valid for 1D-histograms.
596// If the binning is found, but it IsDefault() kTRUE is returned, but
597// no binning is applied.
598//
599Bool_t MH::ApplyBinning(const MParList &plist, TString name, TH1 *h)
600{
601 if (h->InheritsFrom(TH2::Class()) || h->InheritsFrom(TH3::Class()))
602 {
603 gLog << warn << "MH::ApplyBinning: '" << h->GetName() << "' is not a basic TH1 object... no binning applied." << endl;
604 return kFALSE;
605 }
606
607 const MBinning *bins = (MBinning*)plist.FindObject("Binning"+name);
608 if (!bins)
609 {
610 gLog << inf << "Object 'Binning" << name << "' [MBinning] not found... no binning applied." << endl;
611 return kFALSE;
612 }
613
614 if (bins->IsDefault())
615 return kTRUE;
616
617 SetBinning(h, bins);
618 return kTRUE;
619}
620
621void MH::FindGoodLimits(Int_t nbins, Int_t &newbins, Double_t &xmin, Double_t &xmax, Bool_t isInteger)
622{
623#if ROOT_VERSION_CODE > ROOT_VERSION(3,04,01)
624 THLimitsFinder::OptimizeLimits(nbins, newbins, xmin, xmax, isInteger);
625#else
626//*-*-*-*-*-*-*-*-*Find reasonable bin values*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
627//*-* ==========================
628
629 Double_t dx = 0.1*(xmax-xmin);
630 Double_t umin = xmin - dx;
631 Double_t umax = xmax + dx;
632
633 if (umin < 0 && xmin >= 0)
634 umin = 0;
635
636 if (umax > 0 && xmax <= 0)
637 umax = 0;
638
639 Double_t binlow =0;
640 Double_t binhigh =0;
641 Double_t binwidth=0;
642
643 TGaxis::Optimize(umin, umax, nbins, binlow, binhigh, nbins, binwidth, "");
644
645 if (binwidth <= 0 || binwidth > 1.e+39)
646 {
647 xmin = -1;
648 xmax = 1;
649 }
650 else
651 {
652 xmin = binlow;
653 xmax = binhigh;
654 }
655
656 if (isInteger)
657 {
658 Int_t ixmin = (Int_t)xmin;
659 Int_t ixmax = (Int_t)xmax;
660 Double_t dxmin = (Double_t)ixmin;
661 Double_t dxmax = (Double_t)ixmax;
662
663 xmin = xmin<0 && xmin!=dxmin ? dxmin - 1 : dxmin;
664 xmax = xmax>0 && xmax!=dxmax ? dxmax + 1 : dxmax;
665
666 if (xmin>=xmax)
667 xmax = xmin+1;
668
669 Int_t bw = 1 + (Int_t)((xmax-xmin)/nbins);
670
671 nbins = (Int_t)((xmax-xmin)/bw);
672
673 if (xmin+nbins*bw < xmax)
674 {
675 nbins++;
676 xmax = xmin +nbins*bw;
677 }
678 }
679
680 newbins = nbins;
681#endif
682}
683
684// --------------------------------------------------------------------------
685//
686// Returns the lowest entry in a histogram which is greater than gt (eg >0)
687//
688Double_t MH::GetMinimumGT(const TH1 &h, Double_t gt)
689{
690 Double_t min = FLT_MAX;
691
692 const Int_t nx = h.GetXaxis()->GetNbins();
693 const Int_t ny = h.GetYaxis()->GetNbins();
694 const Int_t nz = h.GetZaxis()->GetNbins();
695
696 for (int iz=1; iz<=nz; iz++)
697 for (int iy=1; iy<=ny; iy++)
698 for (int ix=1; ix<=nx; ix++)
699 {
700 const Double_t v = h.GetBinContent(h.GetBin(ix, iy, iz));
701 if (gt<v && v<min)
702 min = v;
703 }
704 return min;
705}
706
707// --------------------------------------------------------------------------
708//
709// Returns the bin center in a logarithmic scale. If the given bin
710// number is <1 it is set to 1. If it is =GetNbins() it is set to
711// GetNbins()
712//
713Double_t MH::GetBinCenterLog(const TAxis &axe, Int_t nbin)
714{
715 if (nbin>axe.GetNbins())
716 nbin = axe.GetNbins();
717
718 if (nbin<1)
719 nbin = 1;
720
721 const Double_t lo = axe.GetBinLowEdge(nbin);
722 const Double_t hi = axe.GetBinUpEdge(nbin);
723
724 const Double_t val = log10(lo) + log10(hi);
725
726 return pow(10, val/2);
727}
728
729// --------------------------------------------------------------------------
730//
731// Draws a copy of the two given histograms. Uses title as the pad title.
732// Also layout the two statistic boxes and a legend.
733//
734void MH::DrawSameCopy(const TH1 &hist1, const TH1 &hist2, const TString title)
735{
736 //
737 // Draw first histogram
738 //
739 TH1 *h1 = ((TH1&)hist1).DrawCopy();
740 gPad->SetBorderMode(0);
741 gPad->Update();
742
743 // FIXME: Also align max/min with set Maximum/Minimum
744 const Double_t maxbin1 = hist1.GetBinContent(hist1.GetMaximumBin());
745 const Double_t maxbin2 = hist2.GetBinContent(hist2.GetMaximumBin());
746 const Double_t minbin1 = hist1.GetBinContent(hist1.GetMinimumBin());
747 const Double_t minbin2 = hist2.GetBinContent(hist2.GetMinimumBin());
748
749 const Double_t max = TMath::Max(maxbin1, maxbin2);
750 const Double_t min = TMath::Min(minbin1, minbin2);
751
752 h1->SetMaximum(max>0?max*1.05:max*0.95);
753 h1->SetMinimum(max>0?min*0.95:min*1.05);
754
755 TPaveText *t = (TPaveText*)gPad->FindObject("title");
756 if (t)
757 {
758 t->SetName((TString)"MHTitle"); // rename object
759 t->Clear(); // clear old lines
760 t->AddText((TString)" "+title+" "); // add the new title
761 t->SetBit(kCanDelete); // make sure object is deleted
762
763 //
764 // FIXME: This is a stupid workaround to hide the redrawn
765 // (see THistPainter::PaintTitle) title
766 //
767 gPad->Modified(); // indicates a change
768 gPad->Update(); // recreates the original title
769 t->Pop(); // bring our title on top
770 }
771
772 //
773 // Rename first statistics box
774 //
775 TPaveStats *s1 = dynamic_cast<TPaveStats*>(gPad->FindObject("stats"));
776 if (!s1)
777 s1 = dynamic_cast<TPaveStats*>(hist1.GetListOfFunctions()->FindObject("stats"));
778 else
779 s1->SetName((TString)"Stat"+hist1.GetTitle());
780
781 if (s1 && s1->GetX2NDC()>0.95)
782 {
783 const Double_t x1 = s1->GetX1NDC()-0.01;
784 s1->SetX1NDC(x1-(s1->GetX2NDC()-s1->GetX1NDC()));
785 s1->SetX2NDC(x1);
786 }
787
788 //
789 // Draw second histogram
790 //
791 TH1 *h2 = ((TH1&)hist2).DrawCopy("sames");
792 gPad->Update();
793
794 //
795 // Draw Legend
796 //
797 TPaveStats *s2 = dynamic_cast<TPaveStats*>(gPad->FindObject("stats"));
798 if (!s2)
799 s2 = dynamic_cast<TPaveStats*>(hist2.GetListOfFunctions()->FindObject("stats"));
800
801 if (s2)
802 {
803 TLegend &l = *new TLegend(s2->GetX1NDC(),
804 s2->GetY1NDC()-0.015-(s2->GetY2NDC()-s2->GetY1NDC())/2,
805 s2->GetX2NDC(),
806 s2->GetY1NDC()-0.01
807 );
808 l.AddEntry(h1, h1->GetTitle());
809 l.AddEntry(h2, h2->GetTitle());
810 l.SetTextSize(s2->GetTextSize());
811 l.SetTextFont(s2->GetTextFont());
812 l.SetBorderSize(s2->GetBorderSize());
813 l.SetBit(kCanDelete);
814 l.Draw();
815 }
816}
817
818// --------------------------------------------------------------------------
819//
820// Draws the two given histograms. Uses title as the pad title.
821// Also layout the two statistic boxes and a legend.
822//
823void MH::DrawSame(TH1 &hist1, TH1 &hist2, const TString title, Bool_t same)
824{
825 //
826 // Draw first histogram
827 //
828 hist1.Draw(same?"same":"");
829 gPad->SetBorderMode(0);
830 gPad->Update();
831
832 if (hist1.GetEntries()>0 && hist2.GetEntries()>0)
833 {
834 const Double_t maxbin1 = hist1.GetBinContent(hist1.GetMaximumBin());
835 const Double_t maxbin2 = hist2.GetBinContent(hist2.GetMaximumBin());
836 const Double_t minbin1 = hist1.GetBinContent(hist1.GetMinimumBin());
837 const Double_t minbin2 = hist2.GetBinContent(hist2.GetMinimumBin());
838
839 const Double_t max = TMath::Max(maxbin1, maxbin2);
840 const Double_t min = TMath::Min(minbin1, minbin2);
841
842 if (max!=min)
843 {
844 hist1.SetMaximum(max>0?max*1.05:max*0.95);
845 hist1.SetMinimum(max>0?min*0.95:min*1.05);
846 }
847 }
848
849 TPaveText *t = (TPaveText*)gPad->FindObject("title");
850 if (t)
851 {
852 t->SetName((TString)"MHTitle"); // rename object
853 t->Clear(); // clear old lines
854 t->AddText((TString)" "+title+" "); // add the new title
855 t->SetBit(kCanDelete); // make sure object is deleted
856
857 //
858 // FIXME: This is a stupid workaround to hide the redrawn
859 // (see THistPainter::PaintTitle) title
860 //
861 gPad->Modified(); // indicates a change
862 gPad->Update(); // recreates the original title
863 t->Pop(); // bring our title on top
864 }
865
866 //
867 // Rename first statistics box
868 //
869 // Where to get the TPaveStats depends on the root version
870 TPaveStats *s1 = dynamic_cast<TPaveStats*>(gPad->FindObject("stats"));
871 if (!s1)
872 s1 = dynamic_cast<TPaveStats*>(hist1.GetListOfFunctions()->FindObject("stats"));
873 else
874 s1->SetName((TString)"Stat"+hist1.GetTitle());
875
876 if (s1 && s1->GetX2NDC()>0.95)
877 {
878 const Double_t x1 = s1->GetX1NDC()-0.01;
879 s1->SetX1NDC(x1-(s1->GetX2NDC()-s1->GetX1NDC()));
880 s1->SetX2NDC(x1);
881 }
882
883 //
884 // Draw second histogram
885 //
886 hist2.Draw("sames");
887 gPad->Update();
888
889 //
890 // Draw Legend
891 //
892 // Where to get the TPaveStats depends on the root version
893 TPaveStats *s2 = dynamic_cast<TPaveStats*>(gPad->FindObject("stats"));
894 if (!s2)
895 s2 = dynamic_cast<TPaveStats*>(hist2.GetListOfFunctions()->FindObject("stats"));
896
897 if (s2)
898 {
899 TLegend &l = *new TLegend(s2->GetX1NDC(),
900 s2->GetY1NDC()-0.015-(s2->GetY2NDC()-s2->GetY1NDC())/2,
901 s2->GetX2NDC(),
902 s2->GetY1NDC()-0.01
903 );
904 l.AddEntry(&hist1, hist1.GetTitle());
905 l.AddEntry(&hist2, hist2.GetTitle());
906 l.SetTextSize(s2->GetTextSize());
907 l.SetTextFont(s2->GetTextFont());
908 l.SetBorderSize(s2->GetBorderSize());
909 l.SetBit(kCanDelete);
910 l.Draw();
911 }
912}
913
914// --------------------------------------------------------------------------
915//
916// If the opt string contains 'nonew' or gPad is not given NULL is returned.
917// Otherwise the present gPad is returned.
918//
919TVirtualPad *MH::GetNewPad(TString &opt)
920{
921 opt.ToLower();
922
923 if (!opt.Contains("nonew"))
924 return NULL;
925
926 opt.ReplaceAll("nonew", "");
927
928 return gPad;
929}
930
931// --------------------------------------------------------------------------
932//
933// Encapsulate the TObject::Clone such, that a cloned TH1 (or derived)
934// object is not added to the current directory, when cloned.
935//
936TObject *MH::Clone(const char *name) const
937{
938 const Bool_t store = TH1::AddDirectoryStatus();
939
940 TH1::AddDirectory(kFALSE);
941 TObject *o = MParContainer::Clone(name);
942 TH1::AddDirectory(store);
943
944 return o;
945}
946
947// --------------------------------------------------------------------------
948//
949// If the opt string contains 'nonew' or gPad is not given a new canvas
950// with size w/h is created. Otherwise the object is cloned and drawn
951// to the present pad. The kCanDelete bit is set for the clone.
952//
953TObject *MH::DrawClone(Option_t *opt, Int_t w, Int_t h) const
954{
955 TString option(opt);
956
957 TVirtualPad *p = GetNewPad(option);
958 if (!p)
959 p = MakeDefCanvas(this, w, h);
960 else
961 if (!option.Contains("same", TString::kIgnoreCase))
962 p->Clear();
963
964 gROOT->SetSelectedPad(NULL);
965
966 TObject *o = MParContainer::DrawClone(option);
967 o->SetBit(kCanDelete);
968 return o;
969}
970
971// --------------------------------------------------------------------------
972//
973// Check whether a class inheriting from MH overwrites the Draw function
974//
975Bool_t MH::OverwritesDraw(TClass *cls) const
976{
977 if (!cls)
978 cls = IsA();
979
980 //
981 // Check whether we reached the base class MTask
982 //
983 if (TString(cls->GetName())=="MH")
984 return kFALSE;
985
986 //
987 // Check whether the class cls overwrites Draw
988 //
989 if (cls->GetMethodAny("Draw"))
990 return kTRUE;
991
992 //
993 // If the class itself doesn't overload it check all it's base classes
994 //
995 TBaseClass *base=NULL;
996 TIter NextBase(cls->GetListOfBases());
997 while ((base=(TBaseClass*)NextBase()))
998 {
999 if (OverwritesDraw(base->GetClassPointer()))
1000 return kTRUE;
1001 }
1002
1003 return kFALSE;
1004}
1005
1006// --------------------------------------------------------------------------
1007//
1008// Cuts the bins containing only zeros at the edges.
1009//
1010// A new number of bins can be defined with nbins != 0
1011// In the case of nbins == 0, no rebinning will take place
1012//
1013// Returns the new (real) number of bins
1014//
1015Int_t MH::StripZeros(TH1 *h, Int_t nbins)
1016{
1017 TAxis &axe = *h->GetXaxis();
1018
1019 const Int_t min1 = axe.GetFirst();
1020 const Int_t max1 = axe.GetLast();
1021 const Int_t range1 = max1-min1;
1022
1023 //
1024 // Check for useless zeros
1025 //
1026 if (range1 == 0)
1027 return 0;
1028
1029 Int_t min2 = 0;
1030 for (int i=min1; i<=max1; i++)
1031 if (h->GetBinContent(i) != 0)
1032 {
1033 min2 = i;
1034 break;
1035 }
1036
1037 //
1038 // If the histogram consists of zeros only
1039 //
1040 if (min2 == max1)
1041 return 0;
1042
1043 Int_t max2 = 0;
1044 for (int i=max1; i>=min2; i--)
1045 if (h->GetBinContent(i) != 0)
1046 {
1047 max2 = i;
1048 break;
1049 }
1050
1051 //
1052 // Appying TAxis->SetRange before ReBin does not work ...
1053 // But this workaround helps quite fine
1054 //
1055 Axis_t min = h->GetBinLowEdge(min2);
1056 Axis_t max = h->GetBinLowEdge(max2)+h->GetBinWidth(max2);
1057
1058 Int_t nbins2 = max2-min2;
1059 //
1060 // Check for rebinning
1061 //
1062 if (nbins > 0)
1063 {
1064 const Int_t ngroup = (Int_t)(nbins2*h->GetNbinsX()/nbins/(max1-min1));
1065 if (ngroup > 1)
1066 {
1067 h->Rebin(ngroup);
1068 nbins2 /= ngroup;
1069 }
1070 }
1071
1072 Int_t newbins = 0;
1073 FindGoodLimits(nbins2, newbins, min, max, kFALSE);
1074 axe.SetRangeUser(min,max);
1075 return axe.GetLast()-axe.GetFirst();
1076}
1077
1078void MH::ProjectionX(TH1D &dest, const TH2 &src, Int_t firstybin, Int_t lastybin)
1079{
1080 //*-*-*-*-*Project a 2-D histogram into a 1-D histogram along X*-*-*-*-*-*-*
1081 //*-* ====================================================
1082 //
1083 // The projection dest is always of the type TH1D.
1084 // The projection is made from the channels along the Y axis
1085 // ranging from firstybin to lastybin included.
1086 // By default, bins 1 to ny are included
1087 // When all bins are included, the number of entries in the projection
1088 // is set to the number of entries of the 2-D histogram, otherwise
1089 // the number of entries is incremented by 1 for all non empty cells.
1090 //
1091 // if Sumw2() was called for dest, the errors are computed.
1092 //
1093 TAxis &axex = *((TH2&)src).GetXaxis();
1094 TAxis &axey = *((TH2&)src).GetYaxis();
1095
1096 const Int_t nx = axex.GetNbins();
1097 const Int_t ny = axey.GetNbins();
1098 if (firstybin < 0)
1099 firstybin = 1;
1100 if (lastybin > ny)
1101 lastybin = ny;
1102
1103 dest.Reset();
1104 SetBinning(&dest, &axex);
1105
1106 // Create the projection histogram
1107 const Bool_t computeErrors = dest.GetSumw2N() ? 1 : 0;
1108
1109 // Fill the projected histogram
1110 for (Int_t binx=0; binx<=nx+1; binx++)
1111 {
1112 Double_t err2 = 0;
1113 for (Int_t biny=firstybin; biny<=lastybin; biny++)
1114 {
1115 const Double_t cont = src.GetCellContent(binx,biny);
1116 const Double_t err1 = src.GetCellError(binx,biny);
1117 err2 += err1*err1;
1118 if (cont)
1119 dest.Fill(axex.GetBinCenter(binx), cont);
1120 }
1121 if (computeErrors)
1122 dest.SetBinError(binx, TMath::Sqrt(err2));
1123 }
1124 if (firstybin <=1 && lastybin >= ny)
1125 dest.SetEntries(src.GetEntries());
1126}
1127
1128void MH::ProjectionY(TH1D &dest, const TH2 &src, Int_t firstxbin, Int_t lastxbin)
1129{
1130 //*-*-*-*-*Project a 2-D histogram into a 1-D histogram along X*-*-*-*-*-*-*
1131 //*-* ====================================================
1132 //
1133 // The projection dest is always of the type TH1D.
1134 // The projection is made from the channels along the Y axis
1135 // ranging from firstybin to lastybin included.
1136 // By default, bins 1 to ny are included
1137 // When all bins are included, the number of entries in the projection
1138 // is set to the number of entries of the 2-D histogram, otherwise
1139 // the number of entries is incremented by 1 for all non empty cells.
1140 //
1141 // if Sumw2() was called for dest, the errors are computed.
1142 //
1143 TAxis &axex = *((TH2&)src).GetXaxis();
1144 TAxis &axey = *((TH2&)src).GetYaxis();
1145
1146 const Int_t nx = axex.GetNbins();
1147 const Int_t ny = axey.GetNbins();
1148 if (firstxbin < 0)
1149 firstxbin = 1;
1150 if (lastxbin > nx)
1151 lastxbin = nx;
1152
1153 dest.Reset();
1154 SetBinning(&dest, &axey);
1155
1156 // Create the projection histogram
1157 const Bool_t computeErrors = dest.GetSumw2N() ? 1 : 0;
1158
1159 // Fill the projected histogram
1160 for (Int_t biny=0; biny<=ny+1; biny++)
1161 {
1162 Double_t err2 = 0;
1163 for (Int_t binx=firstxbin; binx<=lastxbin; binx++)
1164 {
1165 const Double_t cont = src.GetCellContent(binx,biny);
1166 const Double_t err1 = src.GetCellError(binx,biny);
1167 err2 += err1*err1;
1168 if (cont)
1169 dest.Fill(axey.GetBinCenter(biny), cont);
1170 }
1171 if (computeErrors)
1172 dest.SetBinError(biny, TMath::Sqrt(err2));
1173 }
1174 if (firstxbin <=1 && lastxbin >= nx)
1175 dest.SetEntries(src.GetEntries());
1176}
1177
1178// --------------------------------------------------------------------------
1179//
1180// In contradiction to TPad::FindObject this function searches recursively
1181// in a pad for an object. gPad is the default.
1182//
1183TObject *MH::FindObjectInPad(const char *name, TVirtualPad *pad)
1184{
1185 if (!pad)
1186 pad = gPad;
1187
1188 if (!pad)
1189 return NULL;
1190
1191 TObject *o;
1192
1193 TIter Next(pad->GetListOfPrimitives());
1194 while ((o=Next()))
1195 {
1196 if (!strcmp(o->GetName(), name))
1197 return o;
1198
1199 if (o->InheritsFrom("TPad"))
1200 if ((o = FindObjectInPad(name, (TVirtualPad*)o)))
1201 return o;
1202 }
1203 return NULL;
1204}
1205
1206// --------------------------------------------------------------------------
1207//
1208//
1209//
1210TH1I* MH::ProjectArray(const TArrayF &array, Int_t nbins,
1211 const char* name, const char* title)
1212{
1213 const Int_t size = array.GetSize();
1214
1215 TH1I *h1=0;
1216
1217 //check if histogram with identical name exist
1218 TObject *h1obj = gROOT->FindObject(name);
1219 if (h1obj && h1obj->InheritsFrom("TH1I"))
1220 {
1221 h1 = (TH1I*)h1obj;
1222 h1->Reset();
1223 }
1224
1225 Double_t min = size>0 ? array[0] : 0;
1226 Double_t max = size>0 ? array[0] : 1;
1227
1228 // first loop over array to find the min and max
1229 for (Int_t i=1; i<size;i++)
1230 {
1231 max = TMath::Max((Double_t)array[i], max);
1232 min = TMath::Min((Double_t)array[i], min);
1233 }
1234
1235 Int_t newbins = 0;
1236 FindGoodLimits(nbins, newbins, min, max, kFALSE);
1237
1238 if (!h1)
1239 {
1240 h1 = new TH1I(name, title, nbins, min, max);
1241 h1->SetXTitle("");
1242 h1->SetYTitle("Counts");
1243 h1->SetDirectory(NULL);
1244 }
1245
1246 // Second loop to fill the histogram
1247 for (Int_t i=0;i<size;i++)
1248 h1->Fill(array[i]);
1249
1250 return h1;
1251}
1252
1253// --------------------------------------------------------------------------
1254//
1255//
1256//
1257TH1I* MH::ProjectArray(const TArrayD &array, Int_t nbins, const char* name, const char* title)
1258{
1259 const Int_t size = array.GetSize();
1260
1261 Double_t min = size>0 ? array[0] : 0;
1262 Double_t max = size>0 ? array[0] : 1;
1263
1264 TH1I *h1=0;
1265
1266 //check if histogram with identical name exist
1267 TObject *h1obj = gROOT->FindObject(name);
1268 if (h1obj && h1obj->InheritsFrom("TH1I"))
1269 {
1270 h1 = (TH1I*)h1obj;
1271 h1->Reset();
1272 }
1273
1274 // first loop over array to find the min and max
1275 for (Int_t i=1; i<size;i++)
1276 {
1277 max = TMath::Max(array[i], max);
1278 min = TMath::Min(array[i], min);
1279 }
1280
1281 Int_t newbins = 0;
1282 FindGoodLimits(nbins, newbins, min, max, kFALSE);
1283
1284 if (!h1)
1285 {
1286 h1 = new TH1I(name, title, newbins, min, max);
1287 h1->SetXTitle("");
1288 h1->SetYTitle("Counts");
1289 h1->SetDirectory(NULL);
1290 }
1291
1292 // Second loop to fill the histogram
1293 for (Int_t i=0;i<size;i++)
1294 h1->Fill(array[i]);
1295
1296 return h1;
1297}
1298
1299// --------------------------------------------------------------------------
1300//
1301//
1302//
1303TH1I* MH::ProjectArray(const MArrayF &array, Int_t nbins,
1304 const char* name, const char* title)
1305{
1306 return ProjectArray(TArrayF(array.GetSize(),array.GetArray()), nbins, name, title);
1307}
1308
1309// --------------------------------------------------------------------------
1310//
1311//
1312//
1313TH1I* MH::ProjectArray(const MArrayD &array, Int_t nbins, const char* name, const char* title)
1314{
1315 return ProjectArray(TArrayD(array.GetSize(),array.GetArray()), nbins, name, title);
1316}
1317
1318// --------------------------------------------------------------------------
1319//
1320// This is a workaround for rrot-version <5.13/04 to get correct
1321// binomial errors even if weights have been used, do:
1322// h->Divide(h1, h2, 5, 2, "b");
1323// MH::SetBinomialErrors(*h, *h1, *h2, 5, 2);
1324//
1325// see http://root.cern.ch/phpBB2/viewtopic.php?p=14818
1326// see http://savannah.cern.ch/bugs/?20722
1327//
1328void MH::SetBinomialErrors(TH1 &hres, const TH1 &h1, const TH1 &h2, Double_t c1, Double_t c2)
1329{
1330 for (Int_t binx=0; binx<=hres.GetNbinsX()+1; binx++)
1331 {
1332 const Double_t b1 = h1.GetBinContent(binx);
1333 const Double_t b2 = h2.GetBinContent(binx);
1334 const Double_t e1 = h2.GetBinError(binx);
1335 const Double_t e2 = h1.GetBinError(binx);
1336
1337 //const Double_t w = c2*b2 ? (c1*b1)/(c2*b2) : 0;
1338 //const Double_t rc = ((1-2*w)*e1*e1+w*w*e2*e2)/(b2*b2);
1339
1340 const Double_t c = c2==0 ? 1 : c1/c2;
1341 const Double_t u = b2==0 ? 0 : b1/b2;
1342
1343 const Double_t rc = c*((1-2*u)*e1*e1+u*u*e2*e2)/(b2*b2);
1344
1345 hres.SetBinError(binx, TMath::Sqrt(TMath::Abs(rc)));
1346 }
1347}
1348
1349// --------------------------------------------------------------------------
1350//
1351// See MTask::PrintSkipped
1352//
1353void MH::PrintSkipped(UInt_t n, const char *str)
1354{
1355 *fLog << " " << setw(7) << n << " (";
1356 *fLog << Form("%5.1f", 100.*n/GetNumExecutions());
1357 *fLog << "%) Evts skipped: " << str << endl;
1358}
1359
1360// --------------------------------------------------------------------------
1361//
1362// Calls gStyle->SetPalette. Allowed palettes are:
1363// pretty
1364// deepblue: darkblue -> lightblue
1365// lightblue: black -> blue -> white
1366// greyscale: black -> white
1367// glow1: black -> darkred -> orange -> yellow -> white
1368// glow2:
1369// glowsym: lightblue -> blue -> black -> darkred -> orange -> yellow -> white
1370// redish: darkred -> lightred
1371// bluish: darkblue -> lightblue
1372// small1:
1373//
1374// If the palette name contains 'inv' the order of the colors is inverted.
1375//
1376// The second argument determines the number of colors for the palette.
1377// The default is 50. 'pretty' always has 50 colors.
1378//
1379// (Remark: Drawing 3D object like TH2D with surf3 allows a maximum
1380// of 99 colors)
1381//
1382void MH::SetPalette(TString paletteName, Int_t ncol)
1383{
1384 Bool_t found=kFALSE;
1385
1386 paletteName.ToLower();
1387
1388 const Bool_t inverse = paletteName.Contains("inv");
1389
1390 if (paletteName.Contains("pretty"))
1391 {
1392 gStyle->SetPalette(1, 0);
1393 ncol=50;
1394 found=kTRUE;
1395 }
1396
1397 if (paletteName.Contains("deepblue"))
1398 {
1399 Double_t s[5] = { 0.00, 0.34, 0.61, 0.84, 1.00 };
1400 Double_t r[5] = { 0.00, 0.09, 0.18, 0.09, 0.00 };
1401 Double_t g[5] = { 0.01, 0.02, 0.39, 0.68, 0.97 };
1402 Double_t b[5] = { 0.17, 0.39, 0.62, 0.79, 0.97 };
1403 gStyle->CreateGradientColorTable(5, s, r, g, b, ncol);
1404 found=kTRUE;
1405 }
1406
1407 if (paletteName.Contains("lightblue"))
1408 {
1409 Double_t s[5] = { 0.00, 0.34, 0.61, 0.84, 1.00 };
1410 Double_t r[5] = { 0.00, 0.09, 0.18, 0.09, 0.00 };
1411 Double_t g[5] = { 0.00, 0.02, 0.40, 0.70, 1.00 };
1412 Double_t b[5] = { 0.00, 0.27, 0.51, 0.81, 1.00 };
1413 gStyle->CreateGradientColorTable(5, s, r, g, b, ncol);
1414 found=kTRUE;
1415 }
1416
1417 if (paletteName.Contains("greyscale"))
1418 {
1419 double s[2] = {0.00, 1.00};
1420 double r[2] = {0.00, 1.00};
1421 double g[2] = {0.00, 1.00};
1422 double b[2] = {0.00, 1.00};
1423 gStyle->CreateGradientColorTable(2, s, r, g, b, ncol);
1424 found=kTRUE;
1425 }
1426
1427 if (paletteName.Contains("glow1"))
1428 {
1429 double s[5] = {0., 0.10, 0.45, 0.75, 1.00};
1430 double r[5] = {0., 0.35, 0.85, 1.00, 1.00};
1431 double g[5] = {0., 0.10, 0.20, 0.73, 1.00};
1432 double b[5] = {0., 0.03, 0.06, 0.00, 1.00};
1433 gStyle->CreateGradientColorTable(5, s, r, g, b, ncol);
1434 found=kTRUE;
1435 }
1436
1437 if (paletteName.Contains("glow2"))
1438 {
1439 double s[4] = {0.00, 0.50, 0.75, 1.00};
1440 double r[4] = {0.24, 0.67, 1.00, 1.00};
1441 double g[4] = {0.03, 0.04, 0.80, 1.00};
1442 double b[4] = {0.03, 0.04, 0.00, 1.00};
1443 gStyle->CreateGradientColorTable(4, s, r, g, b, ncol);
1444 found=kTRUE;
1445 }
1446
1447 if (paletteName.Contains("glowsym"))
1448 {
1449 double s[8] = {0.00, 0.17, 0.39, 0.50, 0.55, 0.72, 0.88, 1.00};
1450 double r[8] = {0.09, 0.18, 0.09, 0.00, 0.35, 0.85, 1.00, 1.00};
1451 double g[8] = {0.70, 0.40, 0.02, 0.00, 0.10, 0.20, 0.73, 1.00};
1452 double b[8] = {0.81, 0.51, 0.27, 0.00, 0.03, 0.06, 0.00, 1.00};
1453 gStyle->CreateGradientColorTable(8, s, r, g, b, ncol);
1454 found=kTRUE;
1455 }
1456
1457 if (paletteName.Contains("redish"))
1458 {
1459 double s[3] = {0., 0.5, 1.};
1460 double r[3] = {0., 1.0, 1.};
1461 double g[3] = {0., 0.0, 1.};
1462 double b[3] = {0., 0.0, 1.};
1463 gStyle->CreateGradientColorTable(3, s, r, g, b, ncol);
1464 found=kTRUE;
1465 }
1466
1467 if (paletteName.Contains("bluish"))
1468 {
1469 double s[3] = {0., 0.5, 1.};
1470 double r[3] = {0., 0.0, 1.};
1471 double g[3] = {0., 0.0, 1.};
1472 double b[3] = {0., 1.0, 1.};
1473 gStyle->CreateGradientColorTable(3, s, r, g, b, ncol);
1474 found=kTRUE;
1475 }
1476
1477 if (paletteName.Contains("small1"))
1478 {
1479 double s[4] = {0.00, 0.50, 0.95, 1.};
1480 double r[4] = {0.04, 0.28, 0.98, 1.};
1481 double g[4] = {0.28, 0.93, 0.03, 1.};
1482 double b[4] = {0.79, 0.11, 0.03, 1.};
1483 gStyle->CreateGradientColorTable(4, s, r, g, b, ncol);
1484 found=kTRUE;
1485 }
1486 if (paletteName.Contains("pepe"))
1487 {
1488 double s[5] = {0.0, 0.6, 0.7, 0.9, 1.0 };
1489 double r[5] = {0.1, 0.1, 1.0, 1.0, 1.0 };
1490 double g[5] = {0.1, 0.1, 0.0, 1.0, 1.0 };
1491 double b[5] = {0.2, 0.7, 0.0, 0.3, 0.9 };
1492 gStyle->CreateGradientColorTable(5, s, r, g, b, ncol);
1493 found=kTRUE;
1494 }
1495
1496 if (inverse)
1497 {
1498 TArrayI c(ncol);
1499 for (int i=0; i<ncol; i++)
1500 c[ncol-i-1] = gStyle->GetColorPalette(i);
1501 gStyle->SetPalette(ncol, c.GetArray());
1502 }
1503
1504 if (!found)
1505 gLog << warn << "MH::SetPalette: Palette " << paletteName << " unknown... ignored." << endl;
1506}
Note: See TracBrowser for help on using the repository browser.