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

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