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

Last change on this file since 6491 was 5994, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 38.0 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)
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// Applies a given binning to a 1D-histogram
194//
195void MH::SetBinning(TH1 *h, const MBinning *binsx)
196{
197 //
198 // Another strange behaviour: TAxis::Set deletes the axis title!
199 //
200 TAxis &x = *h->GetXaxis();
201
202#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
203 TString xtitle = x.GetTitle();
204#endif
205
206 // All this is reset by TAxis::Set
207 const TAttAxis att(x);
208 const Bool_t tm(x.GetTimeDisplay());
209 const TString tf(x.GetTimeFormat());
210
211 //
212 // This is a necessary workaround if one wants to set
213 // non-equidistant bins after the initialization
214 // TH1D::fNcells must be set correctly.
215 //
216 h->SetBins(binsx->GetNumBins(), 0, 1);
217
218 //
219 // Set the binning of the current histogram to the binning
220 // in one of the two given histograms
221 //
222 x.Set(binsx->GetNumBins(), binsx->GetEdges());
223
224 // All this is reset by TAxis::Set
225 att.Copy(x);
226 x.SetTimeDisplay(tm);
227 x.SetTimeFormat(tf);
228
229#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
230 x.SetTitle(xtitle);
231#endif
232}
233
234// --------------------------------------------------------------------------
235//
236// Applies given binnings to the two axis of a 2D-histogram
237//
238void MH::SetBinning(TH2 *h, const MBinning *binsx, const MBinning *binsy)
239{
240 TAxis &x = *h->GetXaxis();
241 TAxis &y = *h->GetYaxis();
242
243 //
244 // Another strange behaviour: TAxis::Set deletes the axis title!
245 //
246#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
247 TString xtitle = x.GetTitle();
248 TString ytitle = y.GetTitle();
249#endif
250
251 // All this is reset by TAxis::Set
252 const TAttAxis attx(x);
253 const TAttAxis atty(y);
254 const Bool_t tmx(x.GetTimeDisplay());
255 const Bool_t tmy(y.GetTimeDisplay());
256 const TString tfx(x.GetTimeFormat());
257 const TString tfy(y.GetTimeFormat());
258
259 //
260 // This is a necessary workaround if one wants to set
261 // non-equidistant bins after the initialization
262 // TH1D::fNcells must be set correctly.
263 //
264 h->SetBins(binsx->GetNumBins(), 0, 1,
265 binsy->GetNumBins(), 0, 1);
266
267 //
268 // Set the binning of the current histogram to the binning
269 // in one of the two given histograms
270 //
271 x.Set(binsx->GetNumBins(), binsx->GetEdges());
272 y.Set(binsy->GetNumBins(), binsy->GetEdges());
273
274 // All this is reset by TAxis::Set
275 attx.Copy(x);
276 atty.Copy(y);
277 x.SetTimeDisplay(tmx);
278 y.SetTimeDisplay(tmy);
279 x.SetTimeFormat(tfx);
280 y.SetTimeFormat(tfy);
281
282#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
283 x.SetTitle(xtitle);
284 y.SetTitle(ytitle);
285#endif
286}
287
288// --------------------------------------------------------------------------
289//
290// Applies given binnings to the three axis of a 3D-histogram
291//
292void MH::SetBinning(TH3 *h, const MBinning *binsx, const MBinning *binsy, const MBinning *binsz)
293{
294 //
295 // Another strange behaviour: TAxis::Set deletes the axis title!
296 //
297 TAxis &x = *h->GetXaxis();
298 TAxis &y = *h->GetYaxis();
299 TAxis &z = *h->GetZaxis();
300
301#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
302 TString xtitle = x.GetTitle();
303 TString ytitle = y.GetTitle();
304 TString ztitle = z.GetTitle();
305#endif
306
307 // All this is reset by TAxis::Set
308 const TAttAxis attx(x);
309 const TAttAxis atty(y);
310 const TAttAxis attz(z);
311 const Bool_t tmx(x.GetTimeDisplay());
312 const Bool_t tmy(y.GetTimeDisplay());
313 const Bool_t tmz(z.GetTimeDisplay());
314 const TString tfx(x.GetTimeFormat());
315 const TString tfy(y.GetTimeFormat());
316 const TString tfz(z.GetTimeFormat());
317
318 //
319 // This is a necessary workaround if one wants to set
320 // non-equidistant bins after the initialization
321 // TH1D::fNcells must be set correctly.
322 //
323 h->SetBins(binsx->GetNumBins(), 0, 1,
324 binsy->GetNumBins(), 0, 1,
325 binsz->GetNumBins(), 0, 1);
326
327 //
328 // Set the binning of the current histogram to the binning
329 // in one of the two given histograms
330 //
331 x.Set(binsx->GetNumBins(), binsx->GetEdges());
332 y.Set(binsy->GetNumBins(), binsy->GetEdges());
333 z.Set(binsz->GetNumBins(), binsz->GetEdges());
334
335 // All this is reset by TAxis::Set
336 attx.Copy(x);
337 atty.Copy(y);
338 attz.Copy(z);
339 x.SetTimeDisplay(tmx);
340 y.SetTimeDisplay(tmy);
341 z.SetTimeDisplay(tmz);
342 x.SetTimeFormat(tfx);
343 y.SetTimeFormat(tfy);
344 z.SetTimeFormat(tfz);
345
346#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
347 x.SetTitle(xtitle);
348 y.SetTitle(ytitle);
349 z.SetTitle(ztitle);
350#endif
351}
352
353// --------------------------------------------------------------------------
354//
355// Applies given binning (the n+1 edges) to the axis of a 1D-histogram
356//
357void MH::SetBinning(TH1 *h, const TArrayD &binsx)
358{
359 MBinning bx;
360 bx.SetEdges(binsx);
361 SetBinning(h, &bx);
362}
363
364// --------------------------------------------------------------------------
365//
366// Applies given binning (the n+1 edges) to the two axis of a
367// 2D-histogram
368//
369void MH::SetBinning(TH2 *h, const TArrayD &binsx, const TArrayD &binsy)
370{
371 MBinning bx;
372 MBinning by;
373 bx.SetEdges(binsx);
374 by.SetEdges(binsy);
375 SetBinning(h, &bx, &by);
376}
377
378// --------------------------------------------------------------------------
379//
380// Applies given binning (the n+1 edges) to the three axis of a
381// 3D-histogram
382//
383void MH::SetBinning(TH3 *h, const TArrayD &binsx, const TArrayD &binsy, const TArrayD &binsz)
384{
385 MBinning bx;
386 MBinning by;
387 MBinning bz;
388 bx.SetEdges(binsx);
389 by.SetEdges(binsy);
390 bz.SetEdges(binsz);
391 SetBinning(h, &bx, &by, &bz);
392}
393
394// --------------------------------------------------------------------------
395//
396// Applies the binning of a TAxis (eg from a root histogram) to the axis
397// of a 1D-histogram
398//
399void MH::SetBinning(TH1 *h, const TAxis *binsx)
400{
401 const Int_t nx = binsx->GetNbins();
402
403 TArrayD bx(nx+1);
404 for (int i=0; i<nx; i++) bx[i] = binsx->GetBinLowEdge(i+1);
405 bx[nx] = binsx->GetXmax();
406
407 SetBinning(h, bx);
408}
409
410// --------------------------------------------------------------------------
411//
412// Applies the binnings of the TAxis' (eg from a root histogram) to the
413// two axis' of a 2D-histogram
414//
415void MH::SetBinning(TH2 *h, const TAxis *binsx, const TAxis *binsy)
416{
417 const Int_t nx = binsx->GetNbins();
418 const Int_t ny = binsy->GetNbins();
419
420 TArrayD bx(nx+1);
421 TArrayD by(ny+1);
422 for (int i=0; i<nx; i++) bx[i] = binsx->GetBinLowEdge(i+1);
423 for (int i=0; i<ny; i++) by[i] = binsy->GetBinLowEdge(i+1);
424 bx[nx] = binsx->GetXmax();
425 by[ny] = binsy->GetXmax();
426
427 SetBinning(h, bx, by);
428}
429
430// --------------------------------------------------------------------------
431//
432// Applies the binnings of the TAxis' (eg from a root histogram) to the
433// three axis' of a 3D-histogram
434//
435void MH::SetBinning(TH3 *h, const TAxis *binsx, const TAxis *binsy, const TAxis *binsz)
436{
437 const Int_t nx = binsx->GetNbins();
438 const Int_t ny = binsy->GetNbins();
439 const Int_t nz = binsz->GetNbins();
440
441 TArrayD bx(nx+1);
442 TArrayD by(ny+1);
443 TArrayD bz(nz+1);
444 for (int i=0; i<nx; i++) bx[i] = binsx->GetBinLowEdge(i+1);
445 for (int i=0; i<ny; i++) by[i] = binsy->GetBinLowEdge(i+1);
446 for (int i=0; i<nz; i++) bz[i] = binsz->GetBinLowEdge(i+1);
447 bx[nx] = binsx->GetXmax();
448 by[ny] = binsy->GetXmax();
449 bz[nz] = binsz->GetXmax();
450
451 SetBinning(h, bx, by, bz);
452}
453
454// --------------------------------------------------------------------------
455//
456// Applies the binnings of one root-histogram x to another one h
457// Both histograms must be of the same type: TH1, TH2 or TH3
458//
459void MH::SetBinning(TH1 *h, const TH1 *x)
460{
461 if (h->InheritsFrom(TH3::Class()) && x->InheritsFrom(TH3::Class()))
462 {
463 SetBinning((TH3*)h, ((TH1*)x)->GetXaxis(), ((TH1*)x)->GetYaxis(), ((TH1*)x)->GetZaxis());
464 return;
465 }
466 if (h->InheritsFrom(TH3::Class()) || x->InheritsFrom(TH3::Class()))
467 return;
468 if (h->InheritsFrom(TH2::Class()) && x->InheritsFrom(TH2::Class()))
469 {
470 SetBinning((TH2*)h, ((TH1*)x)->GetXaxis(), ((TH1*)x)->GetYaxis());
471 return;
472 }
473 if (h->InheritsFrom(TH2::Class()) || x->InheritsFrom(TH2::Class()))
474 return;
475 if (h->InheritsFrom(TH1::Class()) && x->InheritsFrom(TH1::Class()))
476 {
477 SetBinning(h, ((TH1*)x)->GetXaxis());
478 return;
479 }
480}
481
482void MH::RemoveFirstBin(TH1 &h)
483{
484 if (h.InheritsFrom(TH2::Class()) || h.InheritsFrom(TH3::Class()))
485 return;
486
487 const Bool_t haserr = h.GetSumw2N()>0;
488
489 const Int_t n0 = h.GetNbinsX();
490 if (n0<2)
491 return;
492
493 TArrayD val(n0-1);
494 TArrayD err(haserr ? n0-1 : 0);
495 for (int i=1; i<n0; i++)
496 {
497 val[i-1] = h.GetBinContent(i+1);
498 if (haserr)
499 err[i-1] = h.GetBinError(i+1);
500 }
501
502 MBinning bins;
503 bins.SetEdges(h, 'x');
504 bins.RemoveFirstEdge();
505 bins.Apply(h);
506
507 h.Reset();
508
509 for (int i=1; i<n0; i++)
510 {
511 h.SetBinContent(i, val[i-1]);
512 if (haserr)
513 h.SetBinError(i, err[i-1]);
514 }
515}
516
517// --------------------------------------------------------------------------
518//
519// Multiplies all entries in a TArrayD by a float f
520//
521void MH::ScaleArray(TArrayD &bins, Double_t f)
522{
523 for (int i=0; i<bins.GetSize(); i++)
524 bins[i] *= f;
525}
526
527// --------------------------------------------------------------------------
528//
529// Scales the binning of a TAxis by a float f
530//
531TArrayD MH::ScaleAxis(TAxis &axe, Double_t f)
532{
533 TArrayD arr(axe.GetNbins()+1);
534
535 for (int i=1; i<=axe.GetNbins()+1; i++)
536 arr[i-1] = axe.GetBinLowEdge(i);
537
538 ScaleArray(arr, f);
539
540 return arr;
541}
542
543// --------------------------------------------------------------------------
544//
545// Scales the binning of one, two or three axis of a histogram by a float f
546//
547void MH::ScaleAxis(TH1 *h, Double_t fx, Double_t fy, Double_t fz)
548{
549 if (h->InheritsFrom(TH3::Class()))
550 {
551 SetBinning((TH3*)h,
552 ScaleAxis(*h->GetXaxis(), fx),
553 ScaleAxis(*h->GetYaxis(), fy),
554 ScaleAxis(*h->GetZaxis(), fz));
555 return;
556 }
557
558 if (h->InheritsFrom(TH2::Class()))
559 {
560 SetBinning((TH2*)h,
561 ScaleAxis(*h->GetXaxis(), fx),
562 ScaleAxis(*h->GetYaxis(), fy));
563 return;
564 }
565
566 if (h->InheritsFrom(TH1::Class()))
567 SetBinning(h, ScaleAxis(*h->GetXaxis(), fx));
568}
569
570// --------------------------------------------------------------------------
571//
572// Tries to find a MBinning container with the name "Binning"+name
573// in the given parameter list. If it was found it is applied to the
574// given histogram. This is only valid for 1D-histograms
575//
576Bool_t MH::ApplyBinning(const MParList &plist, TString name, TH1 *h)
577{
578 if (h->InheritsFrom(TH2::Class()) || h->InheritsFrom(TH3::Class()))
579 {
580 gLog << warn << "MH::ApplyBinning: '" << h->GetName() << "' is not a basic TH1 object... no binning applied." << endl;
581 return kFALSE;
582 }
583
584 const MBinning *bins = (MBinning*)plist.FindObject("Binning"+name);
585 if (!bins)
586 {
587 gLog << inf << "Object 'Binning" << name << "' [MBinning] not found... no binning applied." << endl;
588 return kFALSE;
589 }
590
591 SetBinning(h, bins);
592 return kTRUE;
593}
594
595void MH::FindGoodLimits(Int_t nbins, Int_t &newbins, Double_t &xmin, Double_t &xmax, Bool_t isInteger)
596{
597#if ROOT_VERSION_CODE > ROOT_VERSION(3,04,01)
598 THLimitsFinder::OptimizeLimits(nbins, newbins, xmin, xmax, isInteger);
599#else
600//*-*-*-*-*-*-*-*-*Find reasonable bin values*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
601//*-* ==========================
602
603 Double_t dx = 0.1*(xmax-xmin);
604 Double_t umin = xmin - dx;
605 Double_t umax = xmax + dx;
606
607 if (umin < 0 && xmin >= 0)
608 umin = 0;
609
610 if (umax > 0 && xmax <= 0)
611 umax = 0;
612
613 Double_t binlow =0;
614 Double_t binhigh =0;
615 Double_t binwidth=0;
616
617 TGaxis::Optimize(umin, umax, nbins, binlow, binhigh, nbins, binwidth, "");
618
619 if (binwidth <= 0 || binwidth > 1.e+39)
620 {
621 xmin = -1;
622 xmax = 1;
623 }
624 else
625 {
626 xmin = binlow;
627 xmax = binhigh;
628 }
629
630 if (isInteger)
631 {
632 Int_t ixmin = (Int_t)xmin;
633 Int_t ixmax = (Int_t)xmax;
634 Double_t dxmin = (Double_t)ixmin;
635 Double_t dxmax = (Double_t)ixmax;
636
637 xmin = xmin<0 && xmin!=dxmin ? dxmin - 1 : dxmin;
638 xmax = xmax>0 && xmax!=dxmax ? dxmax + 1 : dxmax;
639
640 if (xmin>=xmax)
641 xmax = xmin+1;
642
643 Int_t bw = 1 + (Int_t)((xmax-xmin)/nbins);
644
645 nbins = (Int_t)((xmax-xmin)/bw);
646
647 if (xmin+nbins*bw < xmax)
648 {
649 nbins++;
650 xmax = xmin +nbins*bw;
651 }
652 }
653
654 newbins = nbins;
655#endif
656}
657
658// --------------------------------------------------------------------------
659//
660// Returns the lowest entry in a histogram which is greater than gt (eg >0)
661//
662Double_t MH::GetMinimumGT(const TH1 &h, Double_t gt)
663{
664 Double_t min = FLT_MAX;
665
666 const TAxis &axex = *((TH1&)h).GetXaxis();
667 const TAxis &axey = *((TH1&)h).GetYaxis();
668 const TAxis &axez = *((TH1&)h).GetZaxis();
669
670 for (int iz=1; iz<=axez.GetNbins(); iz++)
671 for (int iy=1; iy<=axey.GetNbins(); iy++)
672 for (int ix=1; ix<=axex.GetNbins(); ix++)
673 {
674 const Double_t v = h.GetBinContent(h.GetBin(ix, iy, iz));
675 if (gt<v && v<min)
676 min = v;
677 }
678 return min;
679}
680
681// --------------------------------------------------------------------------
682//
683// Returns the bin center in a logarithmic scale. If the given bin
684// number is <1 it is set to 1. If it is =GetNbins() it is set to
685// GetNbins()
686//
687Double_t MH::GetBinCenterLog(const TAxis &axe, Int_t nbin)
688{
689 if (nbin>axe.GetNbins())
690 nbin = axe.GetNbins();
691
692 if (nbin<1)
693 nbin = 1;
694
695 const Double_t lo = axe.GetBinLowEdge(nbin);
696 const Double_t hi = axe.GetBinUpEdge(nbin);
697
698 const Double_t val = log10(lo) + log10(hi);
699
700 return pow(10, val/2);
701}
702
703// --------------------------------------------------------------------------
704//
705// Draws a copy of the two given histograms. Uses title as the pad title.
706// Also layout the two statistic boxes and a legend.
707//
708void MH::DrawSameCopy(const TH1 &hist1, const TH1 &hist2, const TString title)
709{
710 //
711 // Draw first histogram
712 //
713 TH1 *h1 = ((TH1&)hist1).DrawCopy();
714 gPad->SetBorderMode(0);
715 gPad->Update();
716
717 // FIXME: Also align max/min with set Maximum/Minimum
718 const Double_t maxbin1 = hist1.GetBinContent(hist1.GetMaximumBin());
719 const Double_t maxbin2 = hist2.GetBinContent(hist2.GetMaximumBin());
720 const Double_t minbin1 = hist1.GetBinContent(hist1.GetMinimumBin());
721 const Double_t minbin2 = hist2.GetBinContent(hist2.GetMinimumBin());
722
723 const Double_t max = TMath::Max(maxbin1, maxbin2);
724 const Double_t min = TMath::Min(minbin1, minbin2);
725
726 h1->SetMaximum(max>0?max*1.05:max*0.95);
727 h1->SetMinimum(max>0?min*0.95:min*1.05);
728
729 TPaveText *t = (TPaveText*)gPad->FindObject("title");
730 if (t)
731 {
732 t->SetName((TString)"MHTitle"); // rename object
733 t->Clear(); // clear old lines
734 t->AddText((TString)" "+title+" "); // add the new title
735 t->SetBit(kCanDelete); // make sure object is deleted
736
737 //
738 // FIXME: This is a stupid workaround to hide the redrawn
739 // (see THistPainter::PaintTitle) title
740 //
741 gPad->Modified(); // indicates a change
742 gPad->Update(); // recreates the original title
743 t->Pop(); // bring our title on top
744 }
745
746 //
747 // Rename first statistics box
748 //
749 TPaveStats *s1 = (TPaveStats*)gPad->FindObject("stats");
750 if (!s1)
751 s1 = (TPaveStats*)hist1.GetListOfFunctions()->FindObject("stats");
752 else
753 s1->SetName((TString)"Stat"+hist1.GetTitle());
754
755 if (s1 && s1->GetX2NDC()>0.95)
756 {
757 const Double_t x1 = s1->GetX1NDC()-0.01;
758 s1->SetX1NDC(x1-(s1->GetX2NDC()-s1->GetX1NDC()));
759 s1->SetX2NDC(x1);
760 }
761
762 //
763 // Draw second histogram
764 //
765 TH1 *h2 = ((TH1&)hist2).DrawCopy("sames");
766 gPad->Update();
767
768 //
769 // Draw Legend
770 //
771 TPaveStats *s2 = (TPaveStats*)gPad->FindObject("stats");
772 if (!s2)
773 s2 = (TPaveStats*)hist2.GetListOfFunctions()->FindObject("stats");
774
775 if (s2)
776 {
777 TLegend &l = *new TLegend(s2->GetX1NDC(),
778 s2->GetY1NDC()-0.015-(s2->GetY2NDC()-s2->GetY1NDC())/2,
779 s2->GetX2NDC(),
780 s2->GetY1NDC()-0.01
781 );
782 l.AddEntry(h1, h1->GetTitle());
783 l.AddEntry(h2, h2->GetTitle());
784 l.SetTextSize(s2->GetTextSize());
785 l.SetTextFont(s2->GetTextFont());
786 l.SetBorderSize(s2->GetBorderSize());
787 l.SetBit(kCanDelete);
788 l.Draw();
789 }
790}
791
792// --------------------------------------------------------------------------
793//
794// Draws the two given histograms. Uses title as the pad title.
795// Also layout the two statistic boxes and a legend.
796//
797void MH::DrawSame(TH1 &hist1, TH1 &hist2, const TString title, Bool_t same)
798{
799 //
800 // Draw first histogram
801 //
802 hist1.Draw(same?"same":"");
803 gPad->SetBorderMode(0);
804 gPad->Update();
805
806 if (hist1.GetEntries()>0 && hist2.GetEntries()>0)
807 {
808 const Double_t maxbin1 = hist1.GetBinContent(hist1.GetMaximumBin());
809 const Double_t maxbin2 = hist2.GetBinContent(hist2.GetMaximumBin());
810 const Double_t minbin1 = hist1.GetBinContent(hist1.GetMinimumBin());
811 const Double_t minbin2 = hist2.GetBinContent(hist2.GetMinimumBin());
812
813 const Double_t max = TMath::Max(maxbin1, maxbin2);
814 const Double_t min = TMath::Min(minbin1, minbin2);
815
816 if (max!=min)
817 {
818 hist1.SetMaximum(max>0?max*1.05:max*0.95);
819 hist1.SetMinimum(max>0?min*0.95:min*1.05);
820 }
821 }
822
823 TPaveText *t = (TPaveText*)gPad->FindObject("title");
824 if (t)
825 {
826 t->SetName((TString)"MHTitle"); // rename object
827 t->Clear(); // clear old lines
828 t->AddText((TString)" "+title+" "); // add the new title
829 t->SetBit(kCanDelete); // make sure object is deleted
830
831 //
832 // FIXME: This is a stupid workaround to hide the redrawn
833 // (see THistPainter::PaintTitle) title
834 //
835 gPad->Modified(); // indicates a change
836 gPad->Update(); // recreates the original title
837 t->Pop(); // bring our title on top
838 }
839
840 //
841 // Rename first statistics box
842 //
843 // Where to get the TPaveStats depends on the root version
844 TPaveStats *s1 = (TPaveStats*)gPad->FindObject("stats");
845 if (!s1)
846 s1 = (TPaveStats*)hist1.GetListOfFunctions()->FindObject("stats");
847 else
848 s1->SetName((TString)"Stat"+hist1.GetTitle());
849
850 if (s1 && s1->GetX2NDC()>0.95)
851 {
852 const Double_t x1 = s1->GetX1NDC()-0.01;
853 s1->SetX1NDC(x1-(s1->GetX2NDC()-s1->GetX1NDC()));
854 s1->SetX2NDC(x1);
855 }
856
857 //
858 // Draw second histogram
859 //
860 hist2.Draw("sames");
861 gPad->Update();
862
863 //
864 // Draw Legend
865 //
866 // Where to get the TPaveStats depends on the root version
867 TPaveStats *s2 = (TPaveStats*)gPad->FindObject("stats");
868 if (!s2)
869 s2 = (TPaveStats*)hist2.GetListOfFunctions()->FindObject("stats");
870
871 if (s2)
872 {
873 TLegend &l = *new TLegend(s2->GetX1NDC(),
874 s2->GetY1NDC()-0.015-(s2->GetY2NDC()-s2->GetY1NDC())/2,
875 s2->GetX2NDC(),
876 s2->GetY1NDC()-0.01
877 );
878 l.AddEntry(&hist1, hist1.GetTitle());
879 l.AddEntry(&hist2, hist2.GetTitle());
880 l.SetTextSize(s2->GetTextSize());
881 l.SetTextFont(s2->GetTextFont());
882 l.SetBorderSize(s2->GetBorderSize());
883 l.SetBit(kCanDelete);
884 l.Draw();
885 }
886}
887
888// --------------------------------------------------------------------------
889//
890// If the opt string contains 'nonew' or gPad is not given NULL is returned.
891// Otherwise the present gPad is returned.
892//
893TVirtualPad *MH::GetNewPad(TString &opt)
894{
895 opt.ToLower();
896
897 if (!opt.Contains("nonew"))
898 return NULL;
899
900 opt.ReplaceAll("nonew", "");
901
902 return gPad;
903}
904
905// --------------------------------------------------------------------------
906//
907// Encapsulate the TObject::Clone such, that a cloned TH1 (or derived)
908// object is not added to the current directory, when cloned.
909//
910TObject *MH::Clone(const char *name) const
911{
912 const Bool_t store = TH1::AddDirectoryStatus();
913
914 TH1::AddDirectory(kFALSE);
915 TObject *o = MParContainer::Clone(name);
916 TH1::AddDirectory(store);
917
918 return o;
919}
920
921// --------------------------------------------------------------------------
922//
923// If the opt string contains 'nonew' or gPad is not given a new canvas
924// with size w/h is created. Otherwise the object is cloned and drawn
925// to the present pad. The kCanDelete bit is set for the clone.
926//
927TObject *MH::DrawClone(Option_t *opt, Int_t w, Int_t h) const
928{
929 TString option(opt);
930
931 TVirtualPad *p = GetNewPad(option);
932 if (!p)
933 p = MakeDefCanvas(this, w, h);
934 else
935 if (!option.Contains("same", TString::kIgnoreCase))
936 p->Clear();
937
938 gROOT->SetSelectedPad(NULL);
939
940 TObject *o = MParContainer::DrawClone(option);
941 o->SetBit(kCanDelete);
942 return o;
943}
944
945// --------------------------------------------------------------------------
946//
947// Check whether a class inheriting from MH overwrites the Draw function
948//
949Bool_t MH::OverwritesDraw(TClass *cls) const
950{
951 if (!cls)
952 cls = IsA();
953
954 //
955 // Check whether we reached the base class MTask
956 //
957 if (TString(cls->GetName())=="MH")
958 return kFALSE;
959
960 //
961 // Check whether the class cls overwrites Draw
962 //
963 if (cls->GetMethodAny("Draw"))
964 return kTRUE;
965
966 //
967 // If the class itself doesn't overload it check all it's base classes
968 //
969 TBaseClass *base=NULL;
970 TIter NextBase(cls->GetListOfBases());
971 while ((base=(TBaseClass*)NextBase()))
972 {
973 if (OverwritesDraw(base->GetClassPointer()))
974 return kTRUE;
975 }
976
977 return kFALSE;
978}
979
980// --------------------------------------------------------------------------
981//
982// Cuts the bins containing only zeros at the edges.
983//
984// A new number of bins can be defined with nbins != 0
985// In the case of nbins == 0, no rebinning will take place
986//
987// Returns the new (real) number of bins
988//
989Int_t MH::StripZeros(TH1 *h, Int_t nbins)
990{
991 TAxis &axe = *h->GetXaxis();
992
993 const Int_t min1 = axe.GetFirst();
994 const Int_t max1 = axe.GetLast();
995 const Int_t range1 = max1-min1;
996
997 //
998 // Check for useless zeros
999 //
1000 if (range1 == 0)
1001 return 0;
1002
1003 Int_t min2 = 0;
1004 for (int i=min1; i<=max1; i++)
1005 if (h->GetBinContent(i) != 0)
1006 {
1007 min2 = i;
1008 break;
1009 }
1010
1011 //
1012 // If the histogram consists of zeros only
1013 //
1014 if (min2 == max1)
1015 return 0;
1016
1017 Int_t max2 = 0;
1018 for (int i=max1; i>=min2; i--)
1019 if (h->GetBinContent(i) != 0)
1020 {
1021 max2 = i;
1022 break;
1023 }
1024
1025 //
1026 // Appying TAxis->SetRange before ReBin does not work ...
1027 // But this workaround helps quite fine
1028 //
1029 Axis_t min = h->GetBinLowEdge(min2);
1030 Axis_t max = h->GetBinLowEdge(max2)+h->GetBinWidth(max2);
1031
1032 Int_t nbins2 = max2-min2;
1033 //
1034 // Check for rebinning
1035 //
1036 if (nbins > 0)
1037 {
1038 const Int_t ngroup = (Int_t)(nbins2*h->GetNbinsX()/nbins/(max1-min1));
1039 if (ngroup > 1)
1040 {
1041 h->Rebin(ngroup);
1042 nbins2 /= ngroup;
1043 }
1044 }
1045
1046 Int_t newbins = 0;
1047 FindGoodLimits(nbins2, newbins, min, max, kFALSE);
1048 axe.SetRangeUser(min,max);
1049 return axe.GetLast()-axe.GetFirst();
1050}
1051
1052void MH::ProjectionX(TH1D &dest, const TH2 &src, Int_t firstybin, Int_t lastybin)
1053{
1054 //*-*-*-*-*Project a 2-D histogram into a 1-D histogram along X*-*-*-*-*-*-*
1055 //*-* ====================================================
1056 //
1057 // The projection dest is always of the type TH1D.
1058 // The projection is made from the channels along the Y axis
1059 // ranging from firstybin to lastybin included.
1060 // By default, bins 1 to ny are included
1061 // When all bins are included, the number of entries in the projection
1062 // is set to the number of entries of the 2-D histogram, otherwise
1063 // the number of entries is incremented by 1 for all non empty cells.
1064 //
1065 // if Sumw2() was called for dest, the errors are computed.
1066 //
1067 TAxis &axex = *((TH2&)src).GetXaxis();
1068 TAxis &axey = *((TH2&)src).GetYaxis();
1069
1070 const Int_t nx = axex.GetNbins();
1071 const Int_t ny = axey.GetNbins();
1072 if (firstybin < 0)
1073 firstybin = 1;
1074 if (lastybin > ny)
1075 lastybin = ny;
1076
1077 dest.Reset();
1078 SetBinning(&dest, &axex);
1079
1080 // Create the projection histogram
1081 const Bool_t computeErrors = dest.GetSumw2N() ? 1 : 0;
1082
1083 // Fill the projected histogram
1084 for (Int_t binx=0; binx<=nx+1; binx++)
1085 {
1086 Double_t err2 = 0;
1087 for (Int_t biny=firstybin; biny<=lastybin; biny++)
1088 {
1089 const Double_t cont = src.GetCellContent(binx,biny);
1090 const Double_t err1 = src.GetCellError(binx,biny);
1091 err2 += err1*err1;
1092 if (cont)
1093 dest.Fill(axex.GetBinCenter(binx), cont);
1094 }
1095 if (computeErrors)
1096 dest.SetBinError(binx, TMath::Sqrt(err2));
1097 }
1098 if (firstybin <=1 && lastybin >= ny)
1099 dest.SetEntries(src.GetEntries());
1100}
1101
1102void MH::ProjectionY(TH1D &dest, const TH2 &src, Int_t firstxbin, Int_t lastxbin)
1103{
1104 //*-*-*-*-*Project a 2-D histogram into a 1-D histogram along X*-*-*-*-*-*-*
1105 //*-* ====================================================
1106 //
1107 // The projection dest is always of the type TH1D.
1108 // The projection is made from the channels along the Y axis
1109 // ranging from firstybin to lastybin included.
1110 // By default, bins 1 to ny are included
1111 // When all bins are included, the number of entries in the projection
1112 // is set to the number of entries of the 2-D histogram, otherwise
1113 // the number of entries is incremented by 1 for all non empty cells.
1114 //
1115 // if Sumw2() was called for dest, the errors are computed.
1116 //
1117 TAxis &axex = *((TH2&)src).GetXaxis();
1118 TAxis &axey = *((TH2&)src).GetYaxis();
1119
1120 const Int_t nx = axex.GetNbins();
1121 const Int_t ny = axey.GetNbins();
1122 if (firstxbin < 0)
1123 firstxbin = 1;
1124 if (lastxbin > nx)
1125 lastxbin = nx;
1126
1127 dest.Reset();
1128 SetBinning(&dest, &axey);
1129
1130 // Create the projection histogram
1131 const Bool_t computeErrors = dest.GetSumw2N() ? 1 : 0;
1132
1133 // Fill the projected histogram
1134 for (Int_t biny=0; biny<=ny+1; biny++)
1135 {
1136 Double_t err2 = 0;
1137 for (Int_t binx=firstxbin; binx<=lastxbin; binx++)
1138 {
1139 const Double_t cont = src.GetCellContent(binx,biny);
1140 const Double_t err1 = src.GetCellError(binx,biny);
1141 err2 += err1*err1;
1142 if (cont)
1143 dest.Fill(axey.GetBinCenter(biny), cont);
1144 }
1145 if (computeErrors)
1146 dest.SetBinError(biny, TMath::Sqrt(err2));
1147 }
1148 if (firstxbin <=1 && lastxbin >= nx)
1149 dest.SetEntries(src.GetEntries());
1150}
1151
1152// --------------------------------------------------------------------------
1153//
1154// In contradiction to TPad::FindObject this function searches recursively
1155// in a pad for an object. gPad is the default.
1156//
1157TObject *MH::FindObjectInPad(const char *name, TVirtualPad *pad)
1158{
1159 if (!pad)
1160 pad = gPad;
1161
1162 if (!pad)
1163 return NULL;
1164
1165 TObject *o;
1166
1167 TIter Next(pad->GetListOfPrimitives());
1168 while ((o=Next()))
1169 {
1170 if (!strcmp(o->GetName(), name))
1171 return o;
1172
1173 if (o->InheritsFrom("TPad"))
1174 if ((o = FindObjectInPad(name, (TVirtualPad*)o)))
1175 return o;
1176 }
1177 return NULL;
1178}
1179
1180// --------------------------------------------------------------------------
1181//
1182// M.Gaug added this withouz Documentation
1183//
1184TH1I* MH::ProjectArray(const TArrayF &array, Int_t nbins,
1185 const char* name, const char* title)
1186{
1187 const Int_t size = array.GetSize();
1188
1189 TH1I *h1=0;
1190
1191 //check if histogram with identical name exist
1192 TObject *h1obj = gROOT->FindObject(name);
1193 if (h1obj && h1obj->InheritsFrom("TH1I"))
1194 {
1195 h1 = (TH1I*)h1obj;
1196 h1->Reset();
1197 }
1198
1199 Double_t min = size>0 ? array[0] : 0;
1200 Double_t max = size>0 ? array[0] : 1;
1201
1202 // first loop over array to find the min and max
1203 for (Int_t i=1; i<size;i++)
1204 {
1205 max = TMath::Max((Double_t)array[i], max);
1206 min = TMath::Min((Double_t)array[i], min);
1207 }
1208
1209 Int_t newbins = 0;
1210 FindGoodLimits(nbins, newbins, min, max, kFALSE);
1211
1212 if (!h1)
1213 {
1214 h1 = new TH1I(name, title, nbins, min, max);
1215 h1->SetXTitle("");
1216 h1->SetYTitle("Counts");
1217 h1->SetDirectory(NULL);
1218 }
1219
1220 // Second loop to fill the histogram
1221 for (Int_t i=0;i<size;i++)
1222 h1->Fill(array[i]);
1223
1224 return h1;
1225}
1226
1227// --------------------------------------------------------------------------
1228//
1229// M.Gaug added this withouz Documentation
1230//
1231TH1I* MH::ProjectArray(const TArrayD &array, Int_t nbins, const char* name, const char* title)
1232{
1233 const Int_t size = array.GetSize();
1234
1235 Double_t min = size>0 ? array[0] : 0;
1236 Double_t max = size>0 ? array[0] : 1;
1237
1238 TH1I *h1=0;
1239
1240 //check if histogram with identical name exist
1241 TObject *h1obj = gROOT->FindObject(name);
1242 if (h1obj && h1obj->InheritsFrom("TH1I"))
1243 {
1244 h1 = (TH1I*)h1obj;
1245 h1->Reset();
1246 }
1247
1248 // first loop over array to find the min and max
1249 for (Int_t i=1; i<size;i++)
1250 {
1251 max = TMath::Max(array[i], max);
1252 min = TMath::Min(array[i], min);
1253 }
1254
1255 Int_t newbins = 0;
1256 FindGoodLimits(nbins, newbins, min, max, kFALSE);
1257
1258 if (!h1)
1259 {
1260 h1 = new TH1I(name, title, newbins, min, max);
1261 h1->SetXTitle("");
1262 h1->SetYTitle("Counts");
1263 h1->SetDirectory(NULL);
1264 }
1265
1266 // Second loop to fill the histogram
1267 for (Int_t i=0;i<size;i++)
1268 h1->Fill(array[i]);
1269
1270 return h1;
1271}
1272
1273// --------------------------------------------------------------------------
1274//
1275// M.Gaug added this withouz Documentation
1276//
1277TH1I* MH::ProjectArray(const MArrayF &array, Int_t nbins,
1278 const char* name, const char* title)
1279{
1280 return ProjectArray(TArrayF(array.GetSize(),array.GetArray()), nbins, name, title);
1281}
1282
1283// --------------------------------------------------------------------------
1284//
1285// M.Gaug added this withouz Documentation
1286//
1287TH1I* MH::ProjectArray(const MArrayD &array, Int_t nbins, const char* name, const char* title)
1288{
1289 return ProjectArray(TArrayD(array.GetSize(),array.GetArray()), nbins, name, title);
1290}
1291
1292// --------------------------------------------------------------------------
1293//
1294// See MTask::PrintSkipped
1295//
1296void MH::PrintSkipped(UInt_t n, const char *str)
1297{
1298 *fLog << " " << setw(7) << n << " (";
1299 *fLog << Form("%5.1f", 100.*n/GetNumExecutions());
1300 *fLog << "%) Evts skipped: " << str << endl;
1301}
Note: See TracBrowser for help on using the repository browser.