source: tags/Mars-V0.9.5.1/mhbase/MH.cc

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