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

Last change on this file since 8813 was 8709, checked in by tbretz, 17 years ago
*** empty log message ***
File size: 48.8 KB
Line 
1/* ======================================================================== *\
2! $Name: not supported by cvs2svn $:$Id: MH.cc,v 1.36 2007-08-25 15:30:24 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-2007
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
621Bool_t MH::ApplyBinning(const MParList &plist, TString x, TString y, TH2 *h)
622{
623 const MBinning *binsx = (MBinning*)plist.FindObject("Binning"+x);
624 if (!binsx)
625 {
626 gLog << inf << "Object 'Binning" << x << "' [MBinning] not found... no binning applied." << endl;
627 return kFALSE;
628 }
629 const MBinning *binsy = (MBinning*)plist.FindObject("Binning"+y);
630 if (!binsy)
631 {
632 gLog << inf << "Object 'Binning" << y << "' [MBinning] not found... no binning applied." << endl;
633 return kFALSE;
634 }
635
636 if (binsx->IsDefault() && binsy->IsDefault())
637 return kTRUE;
638
639 // -------------------------
640 /*
641 MBinning binsx, binsy, binsz;
642 binsx.SetEdges(fHist, 'x');
643 binsy.SetEdges(fHist, 'y');
644 binsz.SetEdges(fHist, 'z');
645 */
646 // -------------------------
647
648 SetBinning(h, binsx, binsy);
649
650 return kTRUE;
651}
652
653Bool_t MH::ApplyBinning(const MParList &plist, TString x, TString y, TString z, TH3 *h)
654{
655 const MBinning *binsx = (MBinning*)plist.FindObject("Binning"+x);
656 if (!binsx)
657 {
658 gLog << inf << "Object 'Binning" << x << "' [MBinning] not found... no binning applied." << endl;
659 return kFALSE;
660 }
661 const MBinning *binsy = (MBinning*)plist.FindObject("Binning"+y);
662 if (!binsy)
663 {
664 gLog << inf << "Object 'Binning" << y << "' [MBinning] not found... no binning applied." << endl;
665 return kFALSE;
666 }
667 const MBinning *binsz = (MBinning*)plist.FindObject("Binning"+z);
668 if (!binsz)
669 {
670 gLog << inf << "Object 'Binning" << z << "' [MBinning] not found... no binning applied." << endl;
671 return kFALSE;
672 }
673
674 if (binsx->IsDefault() && binsy->IsDefault() && binsz->IsDefault())
675 return kTRUE;
676
677 SetBinning(h, binsx, binsy, binsz);
678 return kTRUE;
679}
680
681void MH::FindGoodLimits(Int_t nbins, Int_t &newbins, Double_t &xmin, Double_t &xmax, Bool_t isInteger)
682{
683#if ROOT_VERSION_CODE > ROOT_VERSION(3,04,01)
684 THLimitsFinder::OptimizeLimits(nbins, newbins, xmin, xmax, isInteger);
685#else
686//*-*-*-*-*-*-*-*-*Find reasonable bin values*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
687//*-* ==========================
688
689 Double_t dx = 0.1*(xmax-xmin);
690 Double_t umin = xmin - dx;
691 Double_t umax = xmax + dx;
692
693 if (umin < 0 && xmin >= 0)
694 umin = 0;
695
696 if (umax > 0 && xmax <= 0)
697 umax = 0;
698
699 Double_t binlow =0;
700 Double_t binhigh =0;
701 Double_t binwidth=0;
702
703 TGaxis::Optimize(umin, umax, nbins, binlow, binhigh, nbins, binwidth, "");
704
705 if (binwidth <= 0 || binwidth > 1.e+39)
706 {
707 xmin = -1;
708 xmax = 1;
709 }
710 else
711 {
712 xmin = binlow;
713 xmax = binhigh;
714 }
715
716 if (isInteger)
717 {
718 Int_t ixmin = (Int_t)xmin;
719 Int_t ixmax = (Int_t)xmax;
720 Double_t dxmin = (Double_t)ixmin;
721 Double_t dxmax = (Double_t)ixmax;
722
723 xmin = xmin<0 && xmin!=dxmin ? dxmin - 1 : dxmin;
724 xmax = xmax>0 && xmax!=dxmax ? dxmax + 1 : dxmax;
725
726 if (xmin>=xmax)
727 xmax = xmin+1;
728
729 Int_t bw = 1 + (Int_t)((xmax-xmin)/nbins);
730
731 nbins = (Int_t)((xmax-xmin)/bw);
732
733 if (xmin+nbins*bw < xmax)
734 {
735 nbins++;
736 xmax = xmin +nbins*bw;
737 }
738 }
739
740 newbins = nbins;
741#endif
742}
743
744// --------------------------------------------------------------------------
745//
746// Returns the lowest entry in a histogram which is greater than gt (eg >0)
747//
748Double_t MH::GetMinimumGT(const TH1 &h, Double_t gt)
749{
750 Double_t min = FLT_MAX;
751
752 const Int_t nx = h.GetXaxis()->GetNbins();
753 const Int_t ny = h.GetYaxis()->GetNbins();
754 const Int_t nz = h.GetZaxis()->GetNbins();
755
756 for (int iz=1; iz<=nz; iz++)
757 for (int iy=1; iy<=ny; iy++)
758 for (int ix=1; ix<=nx; ix++)
759 {
760 const Double_t v = h.GetBinContent(h.GetBin(ix, iy, iz));
761 if (gt<v && v<min)
762 min = v;
763 }
764 return min;
765}
766
767// --------------------------------------------------------------------------
768//
769// Returns the bin center in a logarithmic scale. If the given bin
770// number is <1 it is set to 1. If it is =GetNbins() it is set to
771// GetNbins()
772//
773Double_t MH::GetBinCenterLog(const TAxis &axe, Int_t nbin)
774{
775 if (nbin>axe.GetNbins())
776 nbin = axe.GetNbins();
777
778 if (nbin<1)
779 nbin = 1;
780
781 const Double_t lo = axe.GetBinLowEdge(nbin);
782 const Double_t hi = axe.GetBinUpEdge(nbin);
783
784 const Double_t val = log10(lo) + log10(hi);
785
786 return pow(10, val/2);
787}
788
789// --------------------------------------------------------------------------
790//
791// Draws a copy of the two given histograms. Uses title as the pad title.
792// Also layout the two statistic boxes and a legend.
793//
794void MH::DrawSameCopy(const TH1 &hist1, const TH1 &hist2, const TString title)
795{
796 //
797 // Draw first histogram
798 //
799 TH1 *h1 = ((TH1&)hist1).DrawCopy();
800 gPad->SetBorderMode(0);
801 gPad->Update();
802
803 // FIXME: Also align max/min with set Maximum/Minimum
804 const Double_t maxbin1 = hist1.GetBinContent(hist1.GetMaximumBin());
805 const Double_t maxbin2 = hist2.GetBinContent(hist2.GetMaximumBin());
806 const Double_t minbin1 = hist1.GetBinContent(hist1.GetMinimumBin());
807 const Double_t minbin2 = hist2.GetBinContent(hist2.GetMinimumBin());
808
809 const Double_t max = TMath::Max(maxbin1, maxbin2);
810 const Double_t min = TMath::Min(minbin1, minbin2);
811
812 h1->SetMaximum(max>0?max*1.05:max*0.95);
813 h1->SetMinimum(max>0?min*0.95:min*1.05);
814
815 TPaveText *t = (TPaveText*)gPad->FindObject("title");
816 if (t)
817 {
818 t->SetName((TString)"MHTitle"); // rename object
819 t->Clear(); // clear old lines
820 t->AddText((TString)" "+title+" "); // add the new title
821 t->SetBit(kCanDelete); // make sure object is deleted
822
823 //
824 // FIXME: This is a stupid workaround to hide the redrawn
825 // (see THistPainter::PaintTitle) title
826 //
827 gPad->Modified(); // indicates a change
828 gPad->Update(); // recreates the original title
829 t->Pop(); // bring our title on top
830 }
831
832 //
833 // Rename first statistics box
834 //
835 TPaveStats *s1 = dynamic_cast<TPaveStats*>(gPad->FindObject("stats"));
836 if (!s1)
837 s1 = dynamic_cast<TPaveStats*>(hist1.GetListOfFunctions()->FindObject("stats"));
838 else
839 s1->SetName((TString)"Stat"+hist1.GetTitle());
840
841 if (s1 && s1->GetX2NDC()>0.95)
842 {
843 const Double_t x1 = s1->GetX1NDC()-0.01;
844 s1->SetX1NDC(x1-(s1->GetX2NDC()-s1->GetX1NDC()));
845 s1->SetX2NDC(x1);
846 }
847
848 //
849 // Draw second histogram
850 //
851 TH1 *h2 = ((TH1&)hist2).DrawCopy("sames");
852 gPad->Update();
853
854 //
855 // Draw Legend
856 //
857 TPaveStats *s2 = dynamic_cast<TPaveStats*>(gPad->FindObject("stats"));
858 if (!s2)
859 s2 = dynamic_cast<TPaveStats*>(hist2.GetListOfFunctions()->FindObject("stats"));
860
861 if (s2)
862 {
863 TLegend &l = *new TLegend(s2->GetX1NDC(),
864 s2->GetY1NDC()-0.015-(s2->GetY2NDC()-s2->GetY1NDC())/2,
865 s2->GetX2NDC(),
866 s2->GetY1NDC()-0.01
867 );
868 l.AddEntry(h1, h1->GetTitle());
869 l.AddEntry(h2, h2->GetTitle());
870 l.SetTextSize(s2->GetTextSize());
871 l.SetTextFont(s2->GetTextFont());
872 l.SetBorderSize(s2->GetBorderSize());
873 l.SetBit(kCanDelete);
874 l.Draw();
875 }
876}
877
878// --------------------------------------------------------------------------
879//
880// Draws the two given histograms. Uses title as the pad title.
881// Also layout the two statistic boxes and a legend.
882//
883void MH::DrawSame(TH1 &hist1, TH1 &hist2, const TString title, Bool_t same)
884{
885 //
886 // Draw first histogram
887 //
888 hist1.Draw(same?"same":"");
889 gPad->SetBorderMode(0);
890 gPad->Update();
891
892 if (hist1.GetEntries()>0 && hist2.GetEntries()>0)
893 {
894 const Double_t maxbin1 = hist1.GetBinContent(hist1.GetMaximumBin());
895 const Double_t maxbin2 = hist2.GetBinContent(hist2.GetMaximumBin());
896 const Double_t minbin1 = hist1.GetBinContent(hist1.GetMinimumBin());
897 const Double_t minbin2 = hist2.GetBinContent(hist2.GetMinimumBin());
898
899 const Double_t max = TMath::Max(maxbin1, maxbin2);
900 const Double_t min = TMath::Min(minbin1, minbin2);
901
902 if (max!=min)
903 {
904 hist1.SetMaximum(max>0?max*1.05:max*0.95);
905 hist1.SetMinimum(max>0?min*0.95:min*1.05);
906 }
907 }
908
909 TPaveText *t = (TPaveText*)gPad->FindObject("title");
910 if (t)
911 {
912 t->SetName((TString)"MHTitle"); // rename object
913 t->Clear(); // clear old lines
914 t->AddText((TString)" "+title+" "); // add the new title
915 t->SetBit(kCanDelete); // make sure object is deleted
916
917 //
918 // FIXME: This is a stupid workaround to hide the redrawn
919 // (see THistPainter::PaintTitle) title
920 //
921 gPad->Modified(); // indicates a change
922 gPad->Update(); // recreates the original title
923 t->Pop(); // bring our title on top
924 }
925
926 //
927 // Rename first statistics box
928 //
929 // Where to get the TPaveStats depends on the root version
930 TPaveStats *s1 = dynamic_cast<TPaveStats*>(gPad->FindObject("stats"));
931 if (!s1)
932 s1 = dynamic_cast<TPaveStats*>(hist1.GetListOfFunctions()->FindObject("stats"));
933 else
934 s1->SetName((TString)"Stat"+hist1.GetTitle());
935
936 if (s1 && s1->GetX2NDC()>0.95)
937 {
938 const Double_t x1 = s1->GetX1NDC()-0.01;
939 s1->SetX1NDC(x1-(s1->GetX2NDC()-s1->GetX1NDC()));
940 s1->SetX2NDC(x1);
941 }
942
943 //
944 // Draw second histogram
945 //
946 hist2.Draw("sames");
947 gPad->Update();
948
949 //
950 // Draw Legend
951 //
952 // Where to get the TPaveStats depends on the root version
953 TPaveStats *s2 = dynamic_cast<TPaveStats*>(gPad->FindObject("stats"));
954 if (!s2)
955 s2 = dynamic_cast<TPaveStats*>(hist2.GetListOfFunctions()->FindObject("stats"));
956
957 if (s2)
958 {
959 TLegend &l = *new TLegend(s2->GetX1NDC(),
960 s2->GetY1NDC()-0.015-(s2->GetY2NDC()-s2->GetY1NDC())/2,
961 s2->GetX2NDC(),
962 s2->GetY1NDC()-0.01
963 );
964 l.AddEntry(&hist1, hist1.GetTitle());
965 l.AddEntry(&hist2, hist2.GetTitle());
966 l.SetTextSize(s2->GetTextSize());
967 l.SetTextFont(s2->GetTextFont());
968 l.SetBorderSize(s2->GetBorderSize());
969 l.SetBit(kCanDelete);
970 l.Draw();
971 }
972}
973
974// --------------------------------------------------------------------------
975//
976// If the opt string contains 'nonew' or gPad is not given NULL is returned.
977// Otherwise the present gPad is returned.
978//
979TVirtualPad *MH::GetNewPad(TString &opt)
980{
981 opt.ToLower();
982
983 if (!opt.Contains("nonew"))
984 return NULL;
985
986 opt.ReplaceAll("nonew", "");
987
988 return gPad;
989}
990
991// --------------------------------------------------------------------------
992//
993// Encapsulate the TObject::Clone such, that a cloned TH1 (or derived)
994// object is not added to the current directory, when cloned.
995//
996TObject *MH::Clone(const char *name) const
997{
998 const Bool_t store = TH1::AddDirectoryStatus();
999
1000 TH1::AddDirectory(kFALSE);
1001 TObject *o = MParContainer::Clone(name);
1002 TH1::AddDirectory(store);
1003
1004 return o;
1005}
1006
1007// --------------------------------------------------------------------------
1008//
1009// If the opt string contains 'nonew' or gPad is not given a new canvas
1010// with size w/h is created. Otherwise the object is cloned and drawn
1011// to the present pad. The kCanDelete bit is set for the clone.
1012//
1013TObject *MH::DrawClone(Option_t *opt, Int_t w, Int_t h) const
1014{
1015 TString option(opt);
1016
1017 TVirtualPad *p = GetNewPad(option);
1018 if (!p)
1019 p = MakeDefCanvas(this, w, h);
1020 else
1021 if (!option.Contains("same", TString::kIgnoreCase))
1022 p->Clear();
1023
1024 gROOT->SetSelectedPad(NULL);
1025
1026 TObject *o = MParContainer::DrawClone(option);
1027 o->SetBit(kCanDelete);
1028 return o;
1029}
1030
1031// --------------------------------------------------------------------------
1032//
1033// Check whether a class inheriting from MH overwrites the Draw function
1034//
1035Bool_t MH::OverwritesDraw(TClass *cls) const
1036{
1037 if (!cls)
1038 cls = IsA();
1039
1040 //
1041 // Check whether we reached the base class MTask
1042 //
1043 if (TString(cls->GetName())=="MH")
1044 return kFALSE;
1045
1046 //
1047 // Check whether the class cls overwrites Draw
1048 //
1049 if (cls->GetMethodAny("Draw"))
1050 return kTRUE;
1051
1052 //
1053 // If the class itself doesn't overload it check all it's base classes
1054 //
1055 TBaseClass *base=NULL;
1056 TIter NextBase(cls->GetListOfBases());
1057 while ((base=(TBaseClass*)NextBase()))
1058 {
1059 if (OverwritesDraw(base->GetClassPointer()))
1060 return kTRUE;
1061 }
1062
1063 return kFALSE;
1064}
1065
1066// --------------------------------------------------------------------------
1067//
1068// Cuts the bins containing only zeros at the edges.
1069//
1070// A new number of bins can be defined with nbins != 0
1071// In the case of nbins == 0, no rebinning will take place
1072//
1073// Returns the new (real) number of bins
1074//
1075Int_t MH::StripZeros(TH1 *h, Int_t nbins)
1076{
1077 TAxis &axe = *h->GetXaxis();
1078
1079 const Int_t min1 = axe.GetFirst();
1080 const Int_t max1 = axe.GetLast();
1081 const Int_t range1 = max1-min1;
1082
1083 //
1084 // Check for useless zeros
1085 //
1086 if (range1 == 0)
1087 return 0;
1088
1089 Int_t min2 = 0;
1090 for (int i=min1; i<=max1; i++)
1091 if (h->GetBinContent(i) != 0)
1092 {
1093 min2 = i;
1094 break;
1095 }
1096
1097 //
1098 // If the histogram consists of zeros only
1099 //
1100 if (min2 == max1)
1101 return 0;
1102
1103 Int_t max2 = 0;
1104 for (int i=max1; i>=min2; i--)
1105 if (h->GetBinContent(i) != 0)
1106 {
1107 max2 = i;
1108 break;
1109 }
1110
1111 //
1112 // Appying TAxis->SetRange before ReBin does not work ...
1113 // But this workaround helps quite fine
1114 //
1115 Axis_t min = h->GetBinLowEdge(min2);
1116 Axis_t max = h->GetBinLowEdge(max2)+h->GetBinWidth(max2);
1117
1118 Int_t nbins2 = max2-min2;
1119 //
1120 // Check for rebinning
1121 //
1122 if (nbins > 0)
1123 {
1124 const Int_t ngroup = (Int_t)(nbins2*h->GetNbinsX()/nbins/(max1-min1));
1125 if (ngroup > 1)
1126 {
1127 h->Rebin(ngroup);
1128 nbins2 /= ngroup;
1129 }
1130 }
1131
1132 Int_t newbins = 0;
1133 FindGoodLimits(nbins2, newbins, min, max, kFALSE);
1134 axe.SetRangeUser(min,max);
1135 return axe.GetLast()-axe.GetFirst();
1136}
1137
1138void MH::ProjectionX(TH1D &dest, const TH2 &src, Int_t firstybin, Int_t lastybin)
1139{
1140 //*-*-*-*-*Project a 2-D histogram into a 1-D histogram along X*-*-*-*-*-*-*
1141 //*-* ====================================================
1142 //
1143 // The projection dest is always of the type TH1D.
1144 // The projection is made from the channels along the Y axis
1145 // ranging from firstybin to lastybin included.
1146 // By default, bins 1 to ny are included
1147 // When all bins are included, the number of entries in the projection
1148 // is set to the number of entries of the 2-D histogram, otherwise
1149 // the number of entries is incremented by 1 for all non empty cells.
1150 //
1151 // if Sumw2() was called for dest, the errors are computed.
1152 //
1153 TAxis &axex = *((TH2&)src).GetXaxis();
1154 TAxis &axey = *((TH2&)src).GetYaxis();
1155
1156 const Int_t nx = axex.GetNbins();
1157 const Int_t ny = axey.GetNbins();
1158 if (firstybin < 0)
1159 firstybin = 1;
1160 if (lastybin > ny)
1161 lastybin = ny;
1162
1163 dest.Reset();
1164 SetBinning(&dest, &axex);
1165
1166 // Create the projection histogram
1167 const Bool_t computeErrors = dest.GetSumw2N() ? 1 : 0;
1168
1169 // Fill the projected histogram
1170 for (Int_t binx=0; binx<=nx+1; binx++)
1171 {
1172 Double_t err2 = 0;
1173 for (Int_t biny=firstybin; biny<=lastybin; biny++)
1174 {
1175 const Double_t cont = src.GetCellContent(binx,biny);
1176 const Double_t err1 = src.GetCellError(binx,biny);
1177 err2 += err1*err1;
1178 if (cont)
1179 dest.Fill(axex.GetBinCenter(binx), cont);
1180 }
1181 if (computeErrors)
1182 dest.SetBinError(binx, TMath::Sqrt(err2));
1183 }
1184 if (firstybin <=1 && lastybin >= ny)
1185 dest.SetEntries(src.GetEntries());
1186}
1187
1188void MH::ProjectionY(TH1D &dest, const TH2 &src, Int_t firstxbin, Int_t lastxbin)
1189{
1190 //*-*-*-*-*Project a 2-D histogram into a 1-D histogram along X*-*-*-*-*-*-*
1191 //*-* ====================================================
1192 //
1193 // The projection dest is always of the type TH1D.
1194 // The projection is made from the channels along the Y axis
1195 // ranging from firstybin to lastybin included.
1196 // By default, bins 1 to ny are included
1197 // When all bins are included, the number of entries in the projection
1198 // is set to the number of entries of the 2-D histogram, otherwise
1199 // the number of entries is incremented by 1 for all non empty cells.
1200 //
1201 // if Sumw2() was called for dest, the errors are computed.
1202 //
1203 TAxis &axex = *((TH2&)src).GetXaxis();
1204 TAxis &axey = *((TH2&)src).GetYaxis();
1205
1206 const Int_t nx = axex.GetNbins();
1207 const Int_t ny = axey.GetNbins();
1208 if (firstxbin < 0)
1209 firstxbin = 1;
1210 if (lastxbin > nx)
1211 lastxbin = nx;
1212
1213 dest.Reset();
1214 SetBinning(&dest, &axey);
1215
1216 // Create the projection histogram
1217 const Bool_t computeErrors = dest.GetSumw2N() ? 1 : 0;
1218
1219 // Fill the projected histogram
1220 for (Int_t biny=0; biny<=ny+1; biny++)
1221 {
1222 Double_t err2 = 0;
1223 for (Int_t binx=firstxbin; binx<=lastxbin; binx++)
1224 {
1225 const Double_t cont = src.GetCellContent(binx,biny);
1226 const Double_t err1 = src.GetCellError(binx,biny);
1227 err2 += err1*err1;
1228 if (cont)
1229 dest.Fill(axey.GetBinCenter(biny), cont);
1230 }
1231 if (computeErrors)
1232 dest.SetBinError(biny, TMath::Sqrt(err2));
1233 }
1234 if (firstxbin <=1 && lastxbin >= nx)
1235 dest.SetEntries(src.GetEntries());
1236}
1237
1238// --------------------------------------------------------------------------
1239//
1240// In contradiction to TPad::FindObject this function searches recursively
1241// in a pad for an object. gPad is the default.
1242//
1243TObject *MH::FindObjectInPad(const char *name, TVirtualPad *pad)
1244{
1245 if (!pad)
1246 pad = gPad;
1247
1248 if (!pad)
1249 return NULL;
1250
1251 TObject *o;
1252
1253 TIter Next(pad->GetListOfPrimitives());
1254 while ((o=Next()))
1255 {
1256 if (!strcmp(o->GetName(), name))
1257 return o;
1258
1259 if (o->InheritsFrom("TPad"))
1260 if ((o = FindObjectInPad(name, (TVirtualPad*)o)))
1261 return o;
1262 }
1263 return NULL;
1264}
1265
1266// --------------------------------------------------------------------------
1267//
1268//
1269//
1270TH1I* MH::ProjectArray(const TArrayF &array, Int_t nbins,
1271 const char* name, const char* title)
1272{
1273 const Int_t size = array.GetSize();
1274
1275 TH1I *h1=0;
1276
1277 //check if histogram with identical name exist
1278 TObject *h1obj = gROOT->FindObject(name);
1279 if (h1obj && h1obj->InheritsFrom("TH1I"))
1280 {
1281 h1 = (TH1I*)h1obj;
1282 h1->Reset();
1283 }
1284
1285 Double_t min = size>0 ? array[0] : 0;
1286 Double_t max = size>0 ? array[0] : 1;
1287
1288 // first loop over array to find the min and max
1289 for (Int_t i=1; i<size;i++)
1290 {
1291 max = TMath::Max((Double_t)array[i], max);
1292 min = TMath::Min((Double_t)array[i], min);
1293 }
1294
1295 Int_t newbins = 0;
1296 FindGoodLimits(nbins, newbins, min, max, kFALSE);
1297
1298 if (!h1)
1299 {
1300 h1 = new TH1I(name, title, nbins, min, max);
1301 h1->SetXTitle("");
1302 h1->SetYTitle("Counts");
1303 h1->SetDirectory(NULL);
1304 }
1305
1306 // Second loop to fill the histogram
1307 for (Int_t i=0;i<size;i++)
1308 h1->Fill(array[i]);
1309
1310 return h1;
1311}
1312
1313// --------------------------------------------------------------------------
1314//
1315//
1316//
1317TH1I* MH::ProjectArray(const TArrayD &array, Int_t nbins, const char* name, const char* title)
1318{
1319 const Int_t size = array.GetSize();
1320
1321 Double_t min = size>0 ? array[0] : 0;
1322 Double_t max = size>0 ? array[0] : 1;
1323
1324 TH1I *h1=0;
1325
1326 //check if histogram with identical name exist
1327 TObject *h1obj = gROOT->FindObject(name);
1328 if (h1obj && h1obj->InheritsFrom("TH1I"))
1329 {
1330 h1 = (TH1I*)h1obj;
1331 h1->Reset();
1332 }
1333
1334 // first loop over array to find the min and max
1335 for (Int_t i=1; i<size;i++)
1336 {
1337 max = TMath::Max(array[i], max);
1338 min = TMath::Min(array[i], min);
1339 }
1340
1341 Int_t newbins = 0;
1342 FindGoodLimits(nbins, newbins, min, max, kFALSE);
1343
1344 if (!h1)
1345 {
1346 h1 = new TH1I(name, title, newbins, min, max);
1347 h1->SetXTitle("");
1348 h1->SetYTitle("Counts");
1349 h1->SetDirectory(NULL);
1350 }
1351
1352 // Second loop to fill the histogram
1353 for (Int_t i=0;i<size;i++)
1354 h1->Fill(array[i]);
1355
1356 return h1;
1357}
1358
1359// --------------------------------------------------------------------------
1360//
1361//
1362//
1363TH1I* MH::ProjectArray(const MArrayF &array, Int_t nbins,
1364 const char* name, const char* title)
1365{
1366 return ProjectArray(TArrayF(array.GetSize(),array.GetArray()), nbins, name, title);
1367}
1368
1369// --------------------------------------------------------------------------
1370//
1371//
1372//
1373TH1I* MH::ProjectArray(const MArrayD &array, Int_t nbins, const char* name, const char* title)
1374{
1375 return ProjectArray(TArrayD(array.GetSize(),array.GetArray()), nbins, name, title);
1376}
1377
1378// --------------------------------------------------------------------------
1379//
1380// This is a workaround for rrot-version <5.13/04 to get correct
1381// binomial errors even if weights have been used, do:
1382// h->Divide(h1, h2, 5, 2, "b");
1383// MH::SetBinomialErrors(*h, *h1, *h2, 5, 2);
1384//
1385// see http://root.cern.ch/phpBB2/viewtopic.php?p=14818
1386// see http://savannah.cern.ch/bugs/?20722
1387//
1388void MH::SetBinomialErrors(TH1 &hres, const TH1 &h1, const TH1 &h2, Double_t c1, Double_t c2)
1389{
1390 for (Int_t binx=0; binx<=hres.GetNbinsX()+1; binx++)
1391 {
1392 const Double_t b1 = h1.GetBinContent(binx);
1393 const Double_t b2 = h2.GetBinContent(binx);
1394 const Double_t e1 = h1.GetBinError(binx);
1395 const Double_t e2 = h2.GetBinError(binx);
1396
1397 //const Double_t w = c2*b2 ? (c1*b1)/(c2*b2) : 0;
1398 //const Double_t rc = ((1-2*w)*e1*e1+w*w*e2*e2)/(b2*b2);
1399
1400 if (b2==0)
1401 {
1402 hres.SetBinError(binx, 0);
1403 continue;
1404 }
1405
1406 const Double_t c = c2==0 ? 1 : c1/c2;
1407 const Double_t u = b2==0 ? 0 : b1/b2;
1408
1409 const Double_t rc = c*((1-2*u)*e1*e1+u*u*e2*e2)/(b2*b2);
1410
1411 hres.SetBinError(binx, TMath::Sqrt(TMath::Abs(rc)));
1412 }
1413}
1414
1415// --------------------------------------------------------------------------
1416//
1417// Return the first and last bin of the histogram which is not 0
1418//
1419void MH::GetRange(const TH1 &h, Int_t &lo, Int_t &hi)
1420{
1421 lo = 0;
1422 hi = 1;
1423
1424 for (int i=1; i<=h.GetNbinsX(); i++)
1425 {
1426 if (lo==0 && h.GetBinContent(i)>0)
1427 lo = i;
1428
1429 if (h.GetBinContent(i)>0)
1430 hi = i;
1431 }
1432}
1433
1434// --------------------------------------------------------------------------
1435//
1436// Return the lower edge of the first and the upper edge of the last bin
1437// of the histogram which is not 0
1438//
1439void MH::GetRangeUser(const TH1 &h, Axis_t &lo, Axis_t &hi)
1440{
1441 Int_t f, l;
1442 GetRange(h, f, l);
1443
1444 lo = h.GetBinLowEdge(f);
1445 hi = h.GetBinLowEdge(l+1);
1446}
1447
1448// --------------------------------------------------------------------------
1449//
1450// Return the first and last column (x-bin) of the histogram which is not 0.
1451// Therefor a proper projection is produced if the argument is a TH2.
1452//
1453// TH3 are not yet supported
1454//
1455void MH::GetRangeX(const TH1 &hist, Int_t &lo, Int_t &hi)
1456{
1457 if (hist.InheritsFrom(TH3::Class()))
1458 return;
1459
1460 if (hist.InheritsFrom(TH2::Class()))
1461 {
1462 TH1 *pro = static_cast<const TH2&>(hist).ProjectionX();
1463 GetRange(*pro, lo, hi);
1464 delete pro;
1465 return;
1466 }
1467
1468 GetRange(hist, lo, hi);
1469}
1470
1471// --------------------------------------------------------------------------
1472//
1473// Return the first and last row (y-bin) of the histogram which is not 0.
1474// Therefor a proper projection is produced if the argument is a TH2.
1475//
1476// TH3 are not yet supported
1477//
1478void MH::GetRangeY(const TH1 &hist, Int_t &lo, Int_t &hi)
1479{
1480 if (hist.InheritsFrom(TH3::Class()))
1481 return;
1482
1483 if (hist.InheritsFrom(TH2::Class()))
1484 {
1485 TH1 *pro = static_cast<const TH2&>(hist).ProjectionY();
1486 GetRange(*pro, lo, hi);
1487 delete pro;
1488 }
1489}
1490
1491// --------------------------------------------------------------------------
1492//
1493// Return the lower edge of the first and the upper edge of the last bin
1494// of the histogram h returned by GetRangeX
1495//
1496void MH::GetRangeUserX(const TH1 &h, Axis_t &lo, Axis_t &hi)
1497{
1498 Int_t f, l;
1499 GetRangeX(h, f, l);
1500
1501 lo = h.GetXaxis()->GetBinLowEdge(f);
1502 hi = h.GetXaxis()->GetBinLowEdge(l+1);
1503}
1504
1505// --------------------------------------------------------------------------
1506//
1507// Return the lower edge of the first and the upper edge of the last bin
1508// of the histogram h returned by GetRangeY
1509//
1510void MH::GetRangeUserY(const TH1 &h, Axis_t &lo, Axis_t &hi)
1511{
1512 Int_t f, l;
1513 GetRangeY(h, f, l);
1514
1515 lo = h.GetYaxis()->GetBinLowEdge(f);
1516 hi = h.GetYaxis()->GetBinLowEdge(l+1);
1517}
1518
1519// --------------------------------------------------------------------------
1520//
1521// See MTask::PrintSkipped
1522//
1523void MH::PrintSkipped(UInt_t n, const char *str)
1524{
1525 *fLog << " " << setw(7) << n << " (";
1526 *fLog << Form("%5.1f", 100.*n/GetNumExecutions());
1527 *fLog << "%) Evts skipped: " << str << endl;
1528}
1529
1530// --------------------------------------------------------------------------
1531//
1532// Calls gStyle->SetPalette. Allowed palettes are:
1533// pretty
1534// deepblue: darkblue -> lightblue
1535// lightblue: black -> blue -> white
1536// greyscale: black -> white
1537// glow1: black -> darkred -> orange -> yellow -> white
1538// glow2:
1539// glowsym: lightblue -> blue -> black -> darkred -> orange -> yellow -> white
1540// redish: darkred -> lightred
1541// bluish: darkblue -> lightblue
1542// small1:
1543//
1544// If the palette name contains 'inv' the order of the colors is inverted.
1545//
1546// The second argument determines the number of colors for the palette.
1547// The default is 50. 'pretty' always has 50 colors.
1548//
1549// (Remark: Drawing 3D object like TH2D with surf3 allows a maximum
1550// of 99 colors)
1551//
1552void MH::SetPalette(TString paletteName, Int_t ncol)
1553{
1554 Bool_t found=kFALSE;
1555
1556 paletteName.ToLower();
1557
1558 const Bool_t inverse = paletteName.Contains("inv");
1559
1560 if (paletteName.Contains("pretty"))
1561 {
1562 gStyle->SetPalette(1, 0);
1563 ncol=50;
1564 found=kTRUE;
1565 }
1566
1567 if (paletteName.Contains("deepblue"))
1568 {
1569 Double_t s[5] = { 0.00, 0.34, 0.61, 0.84, 1.00 };
1570 Double_t r[5] = { 0.00, 0.09, 0.18, 0.09, 0.00 };
1571 Double_t g[5] = { 0.01, 0.02, 0.39, 0.68, 0.97 };
1572 Double_t b[5] = { 0.17, 0.39, 0.62, 0.79, 0.97 };
1573 gStyle->CreateGradientColorTable(5, s, r, g, b, ncol);
1574 found=kTRUE;
1575 }
1576
1577 if (paletteName.Contains("lightblue"))
1578 {
1579 Double_t s[5] = { 0.00, 0.34, 0.61, 0.84, 1.00 };
1580 Double_t r[5] = { 0.00, 0.09, 0.18, 0.09, 0.00 };
1581 Double_t g[5] = { 0.00, 0.02, 0.40, 0.70, 1.00 };
1582 Double_t b[5] = { 0.00, 0.27, 0.51, 0.81, 1.00 };
1583 gStyle->CreateGradientColorTable(5, s, r, g, b, ncol);
1584 found=kTRUE;
1585 }
1586
1587 if (paletteName.Contains("greyscale"))
1588 {
1589 double s[2] = {0.00, 1.00};
1590 double r[2] = {0.00, 1.00};
1591 double g[2] = {0.00, 1.00};
1592 double b[2] = {0.00, 1.00};
1593 gStyle->CreateGradientColorTable(2, s, r, g, b, ncol);
1594 found=kTRUE;
1595 }
1596
1597 if (paletteName.Contains("glow1"))
1598 {
1599 double s[5] = {0., 0.10, 0.45, 0.75, 1.00};
1600 double r[5] = {0., 0.35, 0.85, 1.00, 1.00};
1601 double g[5] = {0., 0.10, 0.20, 0.73, 1.00};
1602 double b[5] = {0., 0.03, 0.06, 0.00, 1.00};
1603 gStyle->CreateGradientColorTable(5, s, r, g, b, ncol);
1604 found=kTRUE;
1605 }
1606
1607 if (paletteName.Contains("glow2"))
1608 {
1609 double s[4] = {0.00, 0.50, 0.75, 1.00};
1610 double r[4] = {0.24, 0.67, 1.00, 1.00};
1611 double g[4] = {0.03, 0.04, 0.80, 1.00};
1612 double b[4] = {0.03, 0.04, 0.00, 1.00};
1613 gStyle->CreateGradientColorTable(4, s, r, g, b, ncol);
1614 found=kTRUE;
1615 }
1616
1617 if (paletteName.Contains("glowsym"))
1618 {
1619 double s[8] = {0.00, 0.17, 0.39, 0.50, 0.55, 0.72, 0.88, 1.00};
1620 double r[8] = {0.09, 0.18, 0.09, 0.00, 0.35, 0.85, 1.00, 1.00};
1621 double g[8] = {0.70, 0.40, 0.02, 0.00, 0.10, 0.20, 0.73, 1.00};
1622 double b[8] = {0.81, 0.51, 0.27, 0.00, 0.03, 0.06, 0.00, 1.00};
1623 gStyle->CreateGradientColorTable(8, s, r, g, b, ncol);
1624 found=kTRUE;
1625 }
1626
1627 if (paletteName.Contains("redish"))
1628 {
1629 double s[3] = {0., 0.5, 1.};
1630 double r[3] = {0., 1.0, 1.};
1631 double g[3] = {0., 0.0, 1.};
1632 double b[3] = {0., 0.0, 1.};
1633 gStyle->CreateGradientColorTable(3, s, r, g, b, ncol);
1634 found=kTRUE;
1635 }
1636
1637 if (paletteName.Contains("bluish"))
1638 {
1639 double s[3] = {0., 0.5, 1.};
1640 double r[3] = {0., 0.0, 1.};
1641 double g[3] = {0., 0.0, 1.};
1642 double b[3] = {0., 1.0, 1.};
1643 gStyle->CreateGradientColorTable(3, s, r, g, b, ncol);
1644 found=kTRUE;
1645 }
1646
1647 if (paletteName.Contains("small1"))
1648 {
1649 double s[4] = {0.00, 0.50, 0.95, 1.};
1650 double r[4] = {0.04, 0.28, 0.98, 1.};
1651 double g[4] = {0.28, 0.93, 0.03, 1.};
1652 double b[4] = {0.79, 0.11, 0.03, 1.};
1653 gStyle->CreateGradientColorTable(4, s, r, g, b, ncol);
1654 found=kTRUE;
1655 }
1656 if (paletteName.Contains("pepe"))
1657 {
1658 double s[5] = {0.0, 0.6, 0.7, 0.9, 1.0 };
1659 double r[5] = {0.1, 0.1, 1.0, 1.0, 1.0 };
1660 double g[5] = {0.1, 0.1, 0.0, 1.0, 1.0 };
1661 double b[5] = {0.2, 0.7, 0.0, 0.3, 0.9 };
1662 gStyle->CreateGradientColorTable(5, s, r, g, b, ncol);
1663 found=kTRUE;
1664 }
1665
1666 if (inverse)
1667 {
1668 TArrayI c(ncol);
1669 for (int i=0; i<ncol; i++)
1670 c[ncol-i-1] = gStyle->GetColorPalette(i);
1671 gStyle->SetPalette(ncol, c.GetArray());
1672 }
1673
1674 if (!found)
1675 gLog << warn << "MH::SetPalette: Palette " << paletteName << " unknown... ignored." << endl;
1676}
Note: See TracBrowser for help on using the repository browser.