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

Last change on this file since 4966 was 4966, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 39.1 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Thomas Bretz 07/2001 <mailto:tbretz@astro.uni-wuerzburg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2002
21!
22!
23\* ======================================================================== */
24
25//////////////////////////////////////////////////////////////////////////////
26// //
27// MH //
28// //
29// This is a base tasks for mars histograms. It defines a common interface //
30// for filling the histograms with events (MH::Fill) which is used by a //
31// common 'filler' And a SetupFill member function which may be used //
32// by MFillH. The idea is: //
33// 1) If your Histogram can become filled by one single container //
34// (like MHHillas) you overload MH::Fill and it gets called with //
35// a pointer to the container with which it should be filled. //
36// //
37// 2) You histogram needs several containers to get filled. Than you //
38// have to overload MH::SetupFill and get the necessary objects from //
39// the parameter list. Use this objects in Fill to fill your //
40// histogram. //
41// //
42// If you want to create your own histogram class the new class must be //
43// derived from MH (instead of the base MParContainer) and you must //
44// the fill function of MH. This is the function which is called to fill //
45// the histogram(s) by the data of a corresponding parameter container. //
46// //
47// Remark: the static member function (eg MakeDefCanvas) can be called //
48// from everywhere using: MH::MakeDefCanvas(...) //
49// //
50//////////////////////////////////////////////////////////////////////////////
51
52#include "MH.h"
53
54#include <TH1.h>
55#include <TH2.h>
56#include <TH3.h>
57#include <TStyle.h> // TStyle::GetScreenFactor
58#include <TGaxis.h>
59#include <TCanvas.h>
60#include <TLegend.h>
61#include <TPaveStats.h>
62#include <TBaseClass.h>
63#if ROOT_VERSION_CODE > ROOT_VERSION(3,04,01)
64#include <THLimitsFinder.h>
65#endif
66
67#include "MLog.h"
68#include "MLogManip.h"
69
70#include "MParList.h"
71#include "MParContainer.h"
72
73#include "MBinning.h"
74
75#include "MArrayD.h"
76#include "MArrayF.h"
77
78ClassImp(MH);
79
80using namespace std;
81
82// --------------------------------------------------------------------------
83//
84// Default Constructor. It sets name and title only. Typically you won't
85// need to change this.
86//
87MH::MH(const char *name, const char *title)
88{
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
480void MH::RemoveFirstBin(TH1 &h)
481{
482 if (h.InheritsFrom(TH2::Class()) || h.InheritsFrom(TH3::Class()))
483 return;
484
485 const Int_t n0 = h.GetNbinsX();
486 if (n0<2)
487 return;
488
489 TArrayD val(n0-1);
490 TArrayD err(n0-1);
491 for (int i=1; i<n0; i++)
492 {
493 val[i-1] = h.GetBinContent(i+1);
494 err[i-1] = h.GetBinError(i+1);
495 }
496
497 MBinning bins;
498 bins.SetEdges(h, 'x');
499 bins.RemoveFirstEdge();
500 bins.Apply(h);
501
502 h.Reset();
503
504 for (int i=1; i<n0; i++)
505 {
506 h.SetBinContent(i, val[i-1]);
507 h.SetBinError(i, err[i-1]);
508 }
509}
510
511// --------------------------------------------------------------------------
512//
513// Multiplies all entries in a TArrayD by a float f
514//
515void MH::ScaleArray(TArrayD &bins, Double_t f)
516{
517 for (int i=0; i<bins.GetSize(); i++)
518 bins[i] *= f;
519}
520
521// --------------------------------------------------------------------------
522//
523// Scales the binning of a TAxis by a float f
524//
525TArrayD MH::ScaleAxis(TAxis &axe, Double_t f)
526{
527 TArrayD arr(axe.GetNbins()+1);
528
529 for (int i=1; i<=axe.GetNbins()+1; i++)
530 arr[i-1] = axe.GetBinLowEdge(i);
531
532 ScaleArray(arr, f);
533
534 return arr;
535}
536
537// --------------------------------------------------------------------------
538//
539// Scales the binning of one, two or three axis of a histogram by a float f
540//
541void MH::ScaleAxis(TH1 *h, Double_t fx, Double_t fy, Double_t fz)
542{
543 if (h->InheritsFrom(TH3::Class()))
544 {
545 SetBinning((TH3*)h,
546 ScaleAxis(*h->GetXaxis(), fx),
547 ScaleAxis(*h->GetYaxis(), fy),
548 ScaleAxis(*h->GetZaxis(), fz));
549 return;
550 }
551
552 if (h->InheritsFrom(TH2::Class()))
553 {
554 SetBinning((TH2*)h,
555 ScaleAxis(*h->GetXaxis(), fx),
556 ScaleAxis(*h->GetYaxis(), fy));
557 return;
558 }
559
560 if (h->InheritsFrom(TH1::Class()))
561 SetBinning(h, ScaleAxis(*h->GetXaxis(), fx));
562}
563
564// --------------------------------------------------------------------------
565//
566// Tries to find a MBinning container with the name "Binning"+name
567// in the given parameter list. If it was found it is applied to the
568// given histogram. This is only valid for 1D-histograms
569//
570Bool_t MH::ApplyBinning(const MParList &plist, TString name, TH1 *h)
571{
572 if (h->InheritsFrom(TH2::Class()) || h->InheritsFrom(TH3::Class()))
573 {
574 gLog << warn << "MH::ApplyBinning: '" << h->GetName() << "' is not a basic TH1 object... no binning applied." << endl;
575 return kFALSE;
576 }
577
578 const MBinning *bins = (MBinning*)plist.FindObject("Binning"+name);
579 if (!bins)
580 {
581 gLog << inf << "Object 'Binning" << name << "' [MBinning] not found... no binning applied." << endl;
582 return kFALSE;
583 }
584
585 SetBinning(h, bins);
586 return kTRUE;
587}
588
589void MH::FindGoodLimits(Int_t nbins, Int_t &newbins, Double_t &xmin, Double_t &xmax, Bool_t isInteger)
590{
591#if ROOT_VERSION_CODE > ROOT_VERSION(3,04,01)
592 THLimitsFinder::OptimizeLimits(nbins, newbins, xmin, xmax, isInteger);
593#else
594//*-*-*-*-*-*-*-*-*Find reasonable bin values*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
595//*-* ==========================
596
597 Double_t dx = 0.1*(xmax-xmin);
598 Double_t umin = xmin - dx;
599 Double_t umax = xmax + dx;
600
601 if (umin < 0 && xmin >= 0)
602 umin = 0;
603
604 if (umax > 0 && xmax <= 0)
605 umax = 0;
606
607 Double_t binlow =0;
608 Double_t binhigh =0;
609 Double_t binwidth=0;
610
611 TGaxis::Optimize(umin, umax, nbins, binlow, binhigh, nbins, binwidth, "");
612
613 if (binwidth <= 0 || binwidth > 1.e+39)
614 {
615 xmin = -1;
616 xmax = 1;
617 }
618 else
619 {
620 xmin = binlow;
621 xmax = binhigh;
622 }
623
624 if (isInteger)
625 {
626 Int_t ixmin = (Int_t)xmin;
627 Int_t ixmax = (Int_t)xmax;
628 Double_t dxmin = (Double_t)ixmin;
629 Double_t dxmax = (Double_t)ixmax;
630
631 xmin = xmin<0 && xmin!=dxmin ? dxmin - 1 : dxmin;
632 xmax = xmax>0 && xmax!=dxmax ? dxmax + 1 : dxmax;
633
634 if (xmin>=xmax)
635 xmax = xmin+1;
636
637 Int_t bw = 1 + (Int_t)((xmax-xmin)/nbins);
638
639 nbins = (Int_t)((xmax-xmin)/bw);
640
641 if (xmin+nbins*bw < xmax)
642 {
643 nbins++;
644 xmax = xmin +nbins*bw;
645 }
646 }
647
648 newbins = nbins;
649#endif
650}
651
652// --------------------------------------------------------------------------
653//
654// Returns the lowest entry in a histogram which is greater than gt (eg >0)
655//
656Double_t MH::GetMinimumGT(const TH1 &h, Double_t gt)
657{
658 Double_t min = FLT_MAX;
659
660 const TAxis &axex = *((TH1&)h).GetXaxis();
661 const TAxis &axey = *((TH1&)h).GetYaxis();
662 const TAxis &axez = *((TH1&)h).GetZaxis();
663
664 for (int iz=1; iz<=axez.GetNbins(); iz++)
665 for (int iy=1; iy<=axey.GetNbins(); iy++)
666 for (int ix=1; ix<=axex.GetNbins(); ix++)
667 {
668 const Double_t v = h.GetBinContent(h.GetBin(ix, iy, iz));
669 if (gt<v && v<min)
670 min = v;
671 }
672 return min;
673}
674
675// --------------------------------------------------------------------------
676//
677// Returns the bin center in a logarithmic scale. If the given bin
678// number is <1 it is set to 1. If it is =GetNbins() it is set to
679// GetNbins()
680//
681Double_t MH::GetBinCenterLog(const TAxis &axe, Int_t nbin)
682{
683 if (nbin>axe.GetNbins())
684 nbin = axe.GetNbins();
685
686 if (nbin<1)
687 nbin = 1;
688
689 const Double_t lo = axe.GetBinLowEdge(nbin);
690 const Double_t hi = axe.GetBinUpEdge(nbin);
691
692 const Double_t val = log10(lo) + log10(hi);
693
694 return pow(10, val/2);
695}
696
697// --------------------------------------------------------------------------
698//
699// Draws a copy of the two given histograms. Uses title as the pad title.
700// Also layout the two statistic boxes and a legend.
701//
702void MH::DrawSameCopy(const TH1 &hist1, const TH1 &hist2, const TString title)
703{
704 //
705 // Draw first histogram
706 //
707 TH1 *h1 = ((TH1&)hist1).DrawCopy();
708 gPad->SetBorderMode(0);
709 gPad->Update();
710
711 // FIXME: Also align max/min with set Maximum/Minimum
712 const Double_t maxbin1 = hist1.GetBinContent(hist1.GetMaximumBin());
713 const Double_t maxbin2 = hist2.GetBinContent(hist2.GetMaximumBin());
714 const Double_t minbin1 = hist1.GetBinContent(hist1.GetMinimumBin());
715 const Double_t minbin2 = hist2.GetBinContent(hist2.GetMinimumBin());
716
717 const Double_t max = TMath::Max(maxbin1, maxbin2);
718 const Double_t min = TMath::Min(minbin1, minbin2);
719
720 h1->SetMaximum(max>0?max*1.05:max*0.95);
721 h1->SetMinimum(max>0?min*0.95:min*1.05);
722
723 TPaveText *t = (TPaveText*)gPad->FindObject("title");
724 if (t)
725 {
726 t->SetName((TString)"MHTitle"); // rename object
727 t->Clear(); // clear old lines
728 t->AddText((TString)" "+title+" "); // add the new title
729 t->SetBit(kCanDelete); // make sure object is deleted
730
731 //
732 // FIXME: This is a stupid workaround to hide the redrawn
733 // (see THistPainter::PaintTitle) title
734 //
735 gPad->Modified(); // indicates a change
736 gPad->Update(); // recreates the original title
737 t->Pop(); // bring our title on top
738 }
739
740 //
741 // Rename first statistics box
742 //
743 TPaveStats *s1 = (TPaveStats*)gPad->FindObject("stats");
744 if (!s1)
745 s1 = (TPaveStats*)hist1.GetListOfFunctions()->FindObject("stats");
746 else
747 s1->SetName((TString)"Stat"+hist1.GetTitle());
748
749 if (s1 && s1->GetX2NDC()>0.95)
750 {
751 const Double_t x1 = s1->GetX1NDC()-0.01;
752 s1->SetX1NDC(x1-(s1->GetX2NDC()-s1->GetX1NDC()));
753 s1->SetX2NDC(x1);
754 }
755
756 //
757 // Draw second histogram
758 //
759 TH1 *h2 = ((TH1&)hist2).DrawCopy("sames");
760 gPad->Update();
761
762 //
763 // Draw Legend
764 //
765 TPaveStats *s2 = (TPaveStats*)gPad->FindObject("stats");
766 if (!s2)
767 s2 = (TPaveStats*)hist2.GetListOfFunctions()->FindObject("stats");
768
769 if (s2)
770 {
771 TLegend &l = *new TLegend(s2->GetX1NDC(),
772 s2->GetY1NDC()-0.015-(s2->GetY2NDC()-s2->GetY1NDC())/2,
773 s2->GetX2NDC(),
774 s2->GetY1NDC()-0.01
775 );
776 l.AddEntry(h1, h1->GetTitle());
777 l.AddEntry(h2, h2->GetTitle());
778 l.SetTextSize(s2->GetTextSize());
779 l.SetTextFont(s2->GetTextFont());
780 l.SetBorderSize(s2->GetBorderSize());
781 l.SetBit(kCanDelete);
782 l.Draw();
783 }
784}
785
786// --------------------------------------------------------------------------
787//
788// Draws the two given histograms. Uses title as the pad title.
789// Also layout the two statistic boxes and a legend.
790//
791void MH::DrawSame(TH1 &hist1, TH1 &hist2, const TString title)
792{
793 //
794 // Draw first histogram
795 //
796 hist1.Draw();
797 gPad->SetBorderMode(0);
798 gPad->Update();
799
800 if (hist1.GetEntries()>0 && hist2.GetEntries()>0)
801 {
802 const Double_t maxbin1 = hist1.GetBinContent(hist1.GetMaximumBin());
803 const Double_t maxbin2 = hist2.GetBinContent(hist2.GetMaximumBin());
804 const Double_t minbin1 = hist1.GetBinContent(hist1.GetMinimumBin());
805 const Double_t minbin2 = hist2.GetBinContent(hist2.GetMinimumBin());
806
807 const Double_t max = TMath::Max(maxbin1, maxbin2);
808 const Double_t min = TMath::Min(minbin1, minbin2);
809
810 if (max!=min)
811 {
812 hist1.SetMaximum(max>0?max*1.05:max*0.95);
813 hist1.SetMinimum(max>0?min*0.95:min*1.05);
814 }
815 }
816
817 TPaveText *t = (TPaveText*)gPad->FindObject("title");
818 if (t)
819 {
820 t->SetName((TString)"MHTitle"); // rename object
821 t->Clear(); // clear old lines
822 t->AddText((TString)" "+title+" "); // add the new title
823 t->SetBit(kCanDelete); // make sure object is deleted
824
825 //
826 // FIXME: This is a stupid workaround to hide the redrawn
827 // (see THistPainter::PaintTitle) title
828 //
829 gPad->Modified(); // indicates a change
830 gPad->Update(); // recreates the original title
831 t->Pop(); // bring our title on top
832 }
833
834 //
835 // Rename first statistics box
836 //
837 // Where to get the TPaveStats depends on the root version
838 TPaveStats *s1 = (TPaveStats*)gPad->FindObject("stats");
839 if (!s1)
840 s1 = (TPaveStats*)hist1.GetListOfFunctions()->FindObject("stats");
841 else
842 s1->SetName((TString)"Stat"+hist1.GetTitle());
843
844 if (s1 && s1->GetX2NDC()>0.95)
845 {
846 const Double_t x1 = s1->GetX1NDC()-0.01;
847 s1->SetX1NDC(x1-(s1->GetX2NDC()-s1->GetX1NDC()));
848 s1->SetX2NDC(x1);
849 }
850
851 //
852 // Draw second histogram
853 //
854 hist2.Draw("sames");
855 gPad->Update();
856
857 //
858 // Draw Legend
859 //
860 // Where to get the TPaveStats depends on the root version
861 TPaveStats *s2 = (TPaveStats*)gPad->FindObject("stats");
862 if (!s2)
863 s2 = (TPaveStats*)hist2.GetListOfFunctions()->FindObject("stats");
864
865 if (s2)
866 {
867 TLegend &l = *new TLegend(s2->GetX1NDC(),
868 s2->GetY1NDC()-0.015-(s2->GetY2NDC()-s2->GetY1NDC())/2,
869 s2->GetX2NDC(),
870 s2->GetY1NDC()-0.01
871 );
872 l.AddEntry(&hist1, hist1.GetTitle());
873 l.AddEntry(&hist2, hist2.GetTitle());
874 l.SetTextSize(s2->GetTextSize());
875 l.SetTextFont(s2->GetTextFont());
876 l.SetBorderSize(s2->GetBorderSize());
877 l.SetBit(kCanDelete);
878 l.Draw();
879 }
880}
881
882// --------------------------------------------------------------------------
883//
884// If the opt string contains 'nonew' or gPad is not given NULL is returned.
885// Otherwise the present gPad is returned.
886//
887TVirtualPad *MH::GetNewPad(TString &opt)
888{
889 opt.ToLower();
890
891 if (!opt.Contains("nonew"))
892 return NULL;
893
894 opt.ReplaceAll("nonew", "");
895
896 return gPad;
897}
898
899// --------------------------------------------------------------------------
900//
901// Encapsulate the TObject::Clone such, that a cloned TH1 (or derived)
902// object is not added to the current directory, when cloned.
903//
904TObject *MH::Clone(const char *name) const
905{
906 const Bool_t store = TH1::AddDirectoryStatus();
907
908 TH1::AddDirectory(kFALSE);
909 TObject *o = MParContainer::Clone(name);
910 TH1::AddDirectory(store);
911
912 return o;
913}
914
915// --------------------------------------------------------------------------
916//
917// If the opt string contains 'nonew' or gPad is not given a new canvas
918// with size w/h is created. Otherwise the object is cloned and drawn
919// to the present pad. The kCanDelete bit is set for the clone.
920//
921TObject *MH::DrawClone(Option_t *opt, Int_t w, Int_t h) const
922{
923 TString option(opt);
924
925 TVirtualPad *p = GetNewPad(option);
926 if (!p)
927 p = MakeDefCanvas(this, w, h);
928 else
929 if (!option.Contains("same", TString::kIgnoreCase))
930 p->Clear();
931
932 gROOT->SetSelectedPad(NULL);
933
934 TObject *o = MParContainer::DrawClone(option);
935 o->SetBit(kCanDelete);
936 return o;
937}
938
939// --------------------------------------------------------------------------
940//
941// Check whether a class inheriting from MH overwrites the Draw function
942//
943Bool_t MH::OverwritesDraw(TClass *cls) const
944{
945 if (!cls)
946 cls = IsA();
947
948 //
949 // Check whether we reached the base class MTask
950 //
951 if (TString(cls->GetName())=="MH")
952 return kFALSE;
953
954 //
955 // Check whether the class cls overwrites Draw
956 //
957 if (cls->GetMethodAny("Draw"))
958 return kTRUE;
959
960 //
961 // If the class itself doesn't overload it check all it's base classes
962 //
963 TBaseClass *base=NULL;
964 TIter NextBase(cls->GetListOfBases());
965 while ((base=(TBaseClass*)NextBase()))
966 {
967 if (OverwritesDraw(base->GetClassPointer()))
968 return kTRUE;
969 }
970
971 return kFALSE;
972}
973
974// --------------------------------------------------------------------------
975//
976// Cuts the bins containing only zeros at the edges.
977//
978// A new number of bins can be defined with nbins != 0
979// In the case of nbins == 0, no rebinning will take place
980//
981// Returns the new (real) number of bins
982//
983Int_t MH::StripZeros(TH1 *h, Int_t nbins)
984{
985 TAxis &axe = *h->GetXaxis();
986
987 const Int_t min1 = axe.GetFirst();
988 const Int_t max1 = axe.GetLast();
989 const Int_t range1 = max1-min1;
990
991 //
992 // Check for useless zeros
993 //
994 if (range1 == 0)
995 return 0;
996
997 Int_t min2 = 0;
998 for (int i=min1; i<=max1; i++)
999 if (h->GetBinContent(i) != 0)
1000 {
1001 min2 = i;
1002 break;
1003 }
1004
1005 //
1006 // If the histogram consists of zeros only
1007 //
1008 if (min2 == max1)
1009 return 0;
1010
1011 Int_t max2 = 0;
1012 for (int i=max1; i>=min2; i--)
1013 if (h->GetBinContent(i) != 0)
1014 {
1015 max2 = i;
1016 break;
1017 }
1018
1019 //
1020 // Appying TAxis->SetRange before ReBin does not work ...
1021 // But this workaround helps quite fine
1022 //
1023 Axis_t min = h->GetBinLowEdge(min2);
1024 Axis_t max = h->GetBinLowEdge(max2)+h->GetBinWidth(max2);
1025
1026 Int_t nbins2 = max2-min2;
1027 //
1028 // Check for rebinning
1029 //
1030 if (nbins > 0)
1031 {
1032 const Int_t ngroup = (Int_t)(nbins2*h->GetNbinsX()/nbins/(max1-min1));
1033 if (ngroup > 1)
1034 {
1035 h->Rebin(ngroup);
1036 nbins2 /= ngroup;
1037 }
1038 }
1039
1040 Int_t newbins = 0;
1041 FindGoodLimits(nbins2, newbins, min, max, kFALSE);
1042 axe.SetRangeUser(min,max);
1043 return axe.GetLast()-axe.GetFirst();
1044}
1045
1046void MH::ProjectionX(TH1D &dest, const TH2 &src, Int_t firstybin, Int_t lastybin)
1047{
1048 //*-*-*-*-*Project a 2-D histogram into a 1-D histogram along X*-*-*-*-*-*-*
1049 //*-* ====================================================
1050 //
1051 // The projection dest is always of the type TH1D.
1052 // The projection is made from the channels along the Y axis
1053 // ranging from firstybin to lastybin included.
1054 // By default, bins 1 to ny are included
1055 // When all bins are included, the number of entries in the projection
1056 // is set to the number of entries of the 2-D histogram, otherwise
1057 // the number of entries is incremented by 1 for all non empty cells.
1058 //
1059 // if Sumw2() was called for dest, the errors are computed.
1060 //
1061 TAxis &axex = *((TH2&)src).GetXaxis();
1062 TAxis &axey = *((TH2&)src).GetYaxis();
1063
1064 const Int_t nx = axex.GetNbins();
1065 const Int_t ny = axey.GetNbins();
1066 if (firstybin < 0)
1067 firstybin = 1;
1068 if (lastybin > ny)
1069 lastybin = ny;
1070
1071 dest.Reset();
1072 SetBinning(&dest, &axex);
1073
1074 // Create the projection histogram
1075 const Bool_t computeErrors = dest.GetSumw2N() ? 1 : 0;
1076
1077 // Fill the projected histogram
1078 for (Int_t binx=0; binx<=nx+1; binx++)
1079 {
1080 Double_t err2 = 0;
1081 for (Int_t biny=firstybin; biny<=lastybin; biny++)
1082 {
1083 const Double_t cont = src.GetCellContent(binx,biny);
1084 const Double_t err1 = src.GetCellError(binx,biny);
1085 err2 += err1*err1;
1086 if (cont)
1087 dest.Fill(axex.GetBinCenter(binx), cont);
1088 }
1089 if (computeErrors)
1090 dest.SetBinError(binx, TMath::Sqrt(err2));
1091 }
1092 if (firstybin <=1 && lastybin >= ny)
1093 dest.SetEntries(src.GetEntries());
1094}
1095
1096void MH::ProjectionY(TH1D &dest, const TH2 &src, Int_t firstxbin, Int_t lastxbin)
1097{
1098 //*-*-*-*-*Project a 2-D histogram into a 1-D histogram along X*-*-*-*-*-*-*
1099 //*-* ====================================================
1100 //
1101 // The projection dest is always of the type TH1D.
1102 // The projection is made from the channels along the Y axis
1103 // ranging from firstybin to lastybin included.
1104 // By default, bins 1 to ny are included
1105 // When all bins are included, the number of entries in the projection
1106 // is set to the number of entries of the 2-D histogram, otherwise
1107 // the number of entries is incremented by 1 for all non empty cells.
1108 //
1109 // if Sumw2() was called for dest, the errors are computed.
1110 //
1111 TAxis &axex = *((TH2&)src).GetXaxis();
1112 TAxis &axey = *((TH2&)src).GetYaxis();
1113
1114 const Int_t nx = axex.GetNbins();
1115 const Int_t ny = axey.GetNbins();
1116 if (firstxbin < 0)
1117 firstxbin = 1;
1118 if (lastxbin > nx)
1119 lastxbin = nx;
1120
1121 dest.Reset();
1122 SetBinning(&dest, &axey);
1123
1124 // Create the projection histogram
1125 const Bool_t computeErrors = dest.GetSumw2N() ? 1 : 0;
1126
1127 // Fill the projected histogram
1128 for (Int_t biny=0; biny<=ny+1; biny++)
1129 {
1130 Double_t err2 = 0;
1131 for (Int_t binx=firstxbin; binx<=lastxbin; binx++)
1132 {
1133 const Double_t cont = src.GetCellContent(binx,biny);
1134 const Double_t err1 = src.GetCellError(binx,biny);
1135 err2 += err1*err1;
1136 if (cont)
1137 dest.Fill(axey.GetBinCenter(biny), cont);
1138 }
1139 if (computeErrors)
1140 dest.SetBinError(biny, TMath::Sqrt(err2));
1141 }
1142 if (firstxbin <=1 && lastxbin >= nx)
1143 dest.SetEntries(src.GetEntries());
1144}
1145
1146// --------------------------------------------------------------------------
1147//
1148// In contradiction to TPad::FindObject this function searches recursively
1149// in a pad for an object. gPad is the default.
1150//
1151TObject *MH::FindObjectInPad(const char *name, TVirtualPad *pad)
1152{
1153 if (!pad)
1154 pad = gPad;
1155
1156 if (!pad)
1157 return NULL;
1158
1159 TObject *o;
1160
1161 TIter Next(pad->GetListOfPrimitives());
1162 while ((o=Next()))
1163 {
1164 if (!strcmp(o->GetName(), name))
1165 return o;
1166
1167 if (o->InheritsFrom("TPad"))
1168 if ((o = FindObjectInPad(name, (TVirtualPad*)o)))
1169 return o;
1170 }
1171 return NULL;
1172}
1173
1174// --------------------------------------------------------------------------
1175//
1176// M.Gaug added this withouz Documentation
1177//
1178TH1I* MH::ProjectArray(const TArrayF &array, Int_t nbins,
1179 const char* name, const char* title)
1180{
1181 const Int_t size = array.GetSize();
1182
1183 TH1I *h1=0;
1184
1185 //check if histogram with identical name exist
1186 TObject *h1obj = gROOT->FindObject(name);
1187 if (h1obj && h1obj->InheritsFrom("TH1I"))
1188 {
1189 h1 = (TH1I*)h1obj;
1190 h1->Reset();
1191 }
1192
1193 Double_t min = size>0 ? array[0] : 0;
1194 Double_t max = size>0 ? array[0] : 1;
1195
1196 // first loop over array to find the min and max
1197 for (Int_t i=1; i<size;i++)
1198 {
1199 max = TMath::Max((Double_t)array[i], max);
1200 min = TMath::Min((Double_t)array[i], min);
1201 }
1202
1203 Int_t newbins = 0;
1204 FindGoodLimits(nbins, newbins, min, max, kFALSE);
1205
1206 if (!h1)
1207 {
1208 h1 = new TH1I(name, title, nbins, min, max);
1209 h1->SetXTitle("");
1210 h1->SetYTitle("Counts");
1211 h1->SetDirectory(gROOT);
1212 }
1213
1214 // Second loop to fill the histogram
1215 for (Int_t i=0;i<size;i++)
1216 h1->Fill(array[i]);
1217
1218 return h1;
1219}
1220
1221// --------------------------------------------------------------------------
1222//
1223// M.Gaug added this withouz Documentation
1224//
1225TH1I* MH::ProjectArray(const TArrayD &array, Int_t nbins, const char* name, const char* title)
1226{
1227 const Int_t size = array.GetSize();
1228 TH1I *h1=0;
1229
1230 //check if histogram with identical name exist
1231 TObject *h1obj = gROOT->FindObject(name);
1232 if (h1obj && h1obj->InheritsFrom("TH1I"))
1233 {
1234 h1 = (TH1I*)h1obj;
1235 h1->Reset();
1236 }
1237
1238 Double_t min = size>0 ? array[0] : 0;
1239 Double_t max = size>0 ? array[0] : 1;
1240
1241 // first loop over array to find the min and max
1242 for (Int_t i=1; i<size;i++)
1243 {
1244 max = TMath::Max(array[i], max);
1245 min = TMath::Min(array[i], min);
1246 }
1247
1248 Int_t newbins = 0;
1249 FindGoodLimits(nbins, newbins, min, max, kFALSE);
1250
1251 if (!h1)
1252 {
1253 h1 = new TH1I(name, title, newbins, min, max);
1254 h1->SetXTitle("");
1255 h1->SetYTitle("Counts");
1256 h1->SetDirectory(gROOT);
1257 }
1258
1259 // Second loop to fill the histogram
1260 for (Int_t i=0;i<size;i++)
1261 h1->Fill(array[i]);
1262
1263 return h1;
1264}
1265
1266// --------------------------------------------------------------------------
1267//
1268// M.Gaug added this withouz Documentation
1269//
1270TH1I* MH::ProjectArray(const MArrayF &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(gROOT);
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// M.Gaug added this withouz Documentation
1316//
1317TH1I* MH::ProjectArray(const MArrayD &array, Int_t nbins, const char* name, const char* title)
1318{
1319 const Int_t size = array.GetSize();
1320 TH1I *h1=0;
1321
1322 //check if histogram with identical name exist
1323 TObject *h1obj = gROOT->FindObject(name);
1324 if (h1obj && h1obj->InheritsFrom("TH1I"))
1325 {
1326 h1 = (TH1I*)h1obj;
1327 h1->Reset();
1328 }
1329
1330 Double_t min = size>0 ? array[0] : 0;
1331 Double_t max = size>0 ? array[0] : 1;
1332
1333 // first loop over array to find the min and max
1334 for (Int_t i=1; i<size;i++)
1335 {
1336 max = TMath::Max(array[i], max);
1337 min = TMath::Min(array[i], min);
1338 }
1339
1340 Int_t newbins = 0;
1341 FindGoodLimits(nbins, newbins, min, max, kFALSE);
1342
1343 if (!h1)
1344 {
1345 h1 = new TH1I(name, title, newbins, min, max);
1346 h1->SetXTitle("");
1347 h1->SetYTitle("Counts");
1348 h1->SetDirectory(gROOT);
1349 }
1350
1351 // Second loop to fill the histogram
1352 for (Int_t i=0;i<size;i++)
1353 h1->Fill(array[i]);
1354
1355 return h1;
1356}
1357
Note: See TracBrowser for help on using the repository browser.