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

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