source: trunk/MagicSoft/Mars/mhist/MH.cc@ 2032

Last change on this file since 2032 was 2032, checked in by tbretz, 21 years ago
*** empty log message ***
File size: 24.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 <TGaxis.h>
58#include <TCanvas.h>
59#include <TLegend.h>
60#include <TPaveStats.h>
61#include <TBaseClass.h>
62#if ROOT_VERSION_CODE > ROOT_VERSION(3,04,01)
63#include <THLimitsFinder.h>
64#endif
65
66#include "MLog.h"
67#include "MLogManip.h"
68
69#include "MParList.h"
70#include "MParContainer.h"
71
72#include "MBinning.h"
73
74ClassImp(MH);
75
76// --------------------------------------------------------------------------
77//
78// Default Constructor. It sets name and title only. Typically you won't
79// need to change this.
80//
81MH::MH(const char *name, const char *title)
82{
83 //
84 // set the name and title of this object
85 //
86 fName = name ? name : "MH";
87 fTitle = title ? title : "Base class for Mars histograms";
88}
89
90// --------------------------------------------------------------------------
91//
92// If you want to use the automatic filling of your derived class you
93// must overload this function. If it is not overloaded you cannot use
94// FillH with this class. The argument is a pointer to a container
95// in your paremeter list which is specified in the MFillH constructor.
96// If you are not going to use it you should at least add
97// Bool_t MH::Fill(const MParContainer *) { return kTRUE; }
98// to your class definition.
99//
100Bool_t MH::Fill(const MParContainer *par, Double_t w)
101{
102 *fLog << warn << GetDescriptor() << ": Fill not overloaded! Can't be used!" << endl;
103 return kFALSE;
104}
105
106// --------------------------------------------------------------------------
107//
108// This virtual function is ment as a generalized interface to retrieve
109// a pointer to a root histogram from the MH-derived class.
110//
111TH1 *MH::GetHistByName(const TString name)
112{
113 *fLog << warn << GetDescriptor() << ": GetHistByName not overloaded! Can't be used!" << endl;
114 return NULL;
115}
116
117// --------------------------------------------------------------------------
118//
119// This is a function which should replace the creation of default
120// canvases like root does. Because this is inconvinient in some aspects.
121// need to change this.
122// You can specify a name for the default canvas and a title. Also
123// width and height can be given.
124// MakeDefCanvas looks for a canvas with the given name. If now name is
125// given the DefCanvasName of root is used. If no such canvas is existing
126// it is created and returned. If such a canvas already exists a new canvas
127// with a name plus anumber is created (the number is calculated by the
128// number of all existing canvases plus one)
129//
130TCanvas *MH::MakeDefCanvas(TString name, const char *title,
131 const UInt_t w, const UInt_t h)
132{
133 const TList *list = (TList*)gROOT->GetListOfCanvases();
134
135 if (name.IsNull())
136 name = gROOT->GetDefCanvasName();
137
138 if (list->FindObject(name))
139 name += Form(" <%d>", list->GetSize()+1);
140
141 return new TCanvas(name, title, w, h);
142}
143
144// --------------------------------------------------------------------------
145//
146// This function works like MakeDefCanvas(name, title, w, h) but name
147// and title are retrieved from the given TObject.
148//
149TCanvas *MH::MakeDefCanvas(const TObject *obj,
150 const UInt_t w, const UInt_t h)
151{
152 return MakeDefCanvas(obj->GetName(), obj->GetTitle(), w, h);
153}
154
155// --------------------------------------------------------------------------
156//
157// Applies a given binning to a 1D-histogram
158//
159void MH::SetBinning(TH1 *h, const MBinning *binsx)
160{
161 //
162 // Another strange behaviour: TAxis::Set deletes the axis title!
163 //
164 TAxis &x = *h->GetXaxis();
165
166#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
167 TString xtitle = x.GetTitle();
168#endif
169
170 //
171 // This is a necessary workaround if one wants to set
172 // non-equidistant bins after the initialization
173 // TH1D::fNcells must be set correctly.
174 //
175 h->SetBins(binsx->GetNumBins(), 0, 1);
176
177 //
178 // Set the binning of the current histogram to the binning
179 // in one of the two given histograms
180 //
181 x.Set(binsx->GetNumBins(), binsx->GetEdges());
182#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
183 x.SetTitle(xtitle);
184#endif
185}
186
187// --------------------------------------------------------------------------
188//
189// Applies given binnings to the two axis of a 2D-histogram
190//
191void MH::SetBinning(TH2 *h, const MBinning *binsx, const MBinning *binsy)
192{
193 TAxis &x = *h->GetXaxis();
194 TAxis &y = *h->GetYaxis();
195
196 //
197 // Another strange behaviour: TAxis::Set deletes the axis title!
198 //
199#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
200 TString xtitle = x.GetTitle();
201 TString ytitle = y.GetTitle();
202#endif
203
204 //
205 // This is a necessary workaround if one wants to set
206 // non-equidistant bins after the initialization
207 // TH1D::fNcells must be set correctly.
208 //
209 h->SetBins(binsx->GetNumBins(), 0, 1,
210 binsy->GetNumBins(), 0, 1);
211
212 //
213 // Set the binning of the current histogram to the binning
214 // in one of the two given histograms
215 //
216 x.Set(binsx->GetNumBins(), binsx->GetEdges());
217 y.Set(binsy->GetNumBins(), binsy->GetEdges());
218#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
219 x.SetTitle(xtitle);
220 y.SetTitle(ytitle);
221#endif
222}
223
224// --------------------------------------------------------------------------
225//
226// Applies given binnings to the three axis of a 3D-histogram
227//
228void MH::SetBinning(TH3 *h, const MBinning *binsx, const MBinning *binsy, const MBinning *binsz)
229{
230 //
231 // Another strange behaviour: TAxis::Set deletes the axis title!
232 //
233 TAxis &x = *h->GetXaxis();
234 TAxis &y = *h->GetYaxis();
235 TAxis &z = *h->GetZaxis();
236
237#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
238 TString xtitle = x.GetTitle();
239 TString ytitle = y.GetTitle();
240 TString ztitle = z.GetTitle();
241#endif
242
243 //
244 // This is a necessary workaround if one wants to set
245 // non-equidistant bins after the initialization
246 // TH1D::fNcells must be set correctly.
247 //
248 h->SetBins(binsx->GetNumBins(), 0, 1,
249 binsy->GetNumBins(), 0, 1,
250 binsz->GetNumBins(), 0, 1);
251
252 //
253 // Set the binning of the current histogram to the binning
254 // in one of the two given histograms
255 //
256 x.Set(binsx->GetNumBins(), binsx->GetEdges());
257 y.Set(binsy->GetNumBins(), binsy->GetEdges());
258 z.Set(binsz->GetNumBins(), binsz->GetEdges());
259#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
260 x.SetTitle(xtitle);
261 y.SetTitle(ytitle);
262 z.SetTitle(ztitle);
263#endif
264}
265
266// --------------------------------------------------------------------------
267//
268// Applies given binning (the n+1 edges) to the axis of a 1D-histogram
269//
270void MH::SetBinning(TH1 *h, const TArrayD &binsx)
271{
272 MBinning bx;
273 bx.SetEdges(binsx);
274 SetBinning(h, &bx);
275}
276
277// --------------------------------------------------------------------------
278//
279// Applies given binning (the n+1 edges) to the two axis of a
280// 2D-histogram
281//
282void MH::SetBinning(TH2 *h, const TArrayD &binsx, const TArrayD &binsy)
283{
284 MBinning bx;
285 MBinning by;
286 bx.SetEdges(binsx);
287 by.SetEdges(binsy);
288 SetBinning(h, &bx, &by);
289}
290
291// --------------------------------------------------------------------------
292//
293// Applies given binning (the n+1 edges) to the three axis of a
294// 3D-histogram
295//
296void MH::SetBinning(TH3 *h, const TArrayD &binsx, const TArrayD &binsy, const TArrayD &binsz)
297{
298 MBinning bx;
299 MBinning by;
300 MBinning bz;
301 bx.SetEdges(binsx);
302 by.SetEdges(binsy);
303 bz.SetEdges(binsz);
304 SetBinning(h, &bx, &by, &bz);
305}
306
307// --------------------------------------------------------------------------
308//
309// Applies the binning of a TAxis (eg from a root histogram) to the axis
310// of a 1D-histogram
311//
312void MH::SetBinning(TH1 *h, const TAxis *binsx)
313{
314 const Int_t nx = binsx->GetNbins();
315
316 TArrayD bx(nx+1);
317 for (int i=0; i<nx; i++) bx[i] = binsx->GetBinLowEdge(i+1);
318 bx[nx] = binsx->GetXmax();
319
320 SetBinning(h, bx);
321}
322
323// --------------------------------------------------------------------------
324//
325// Applies the binnings of the TAxis' (eg from a root histogram) to the
326// two axis' of a 2D-histogram
327//
328void MH::SetBinning(TH2 *h, const TAxis *binsx, const TAxis *binsy)
329{
330 const Int_t nx = binsx->GetNbins();
331 const Int_t ny = binsy->GetNbins();
332
333 TArrayD bx(nx+1);
334 TArrayD by(ny+1);
335 for (int i=0; i<nx; i++) bx[i] = binsx->GetBinLowEdge(i+1);
336 for (int i=0; i<ny; i++) by[i] = binsy->GetBinLowEdge(i+1);
337 bx[nx] = binsx->GetXmax();
338 by[ny] = binsy->GetXmax();
339
340 SetBinning(h, bx, by);
341}
342
343// --------------------------------------------------------------------------
344//
345// Applies the binnings of the TAxis' (eg from a root histogram) to the
346// three axis' of a 3D-histogram
347//
348void MH::SetBinning(TH3 *h, const TAxis *binsx, const TAxis *binsy, const TAxis *binsz)
349{
350 const Int_t nx = binsx->GetNbins();
351 const Int_t ny = binsy->GetNbins();
352 const Int_t nz = binsz->GetNbins();
353
354 TArrayD bx(nx+1);
355 TArrayD by(ny+1);
356 TArrayD bz(nz+1);
357 for (int i=0; i<nx; i++) bx[i] = binsx->GetBinLowEdge(i+1);
358 for (int i=0; i<ny; i++) by[i] = binsy->GetBinLowEdge(i+1);
359 for (int i=0; i<nz; i++) bz[i] = binsz->GetBinLowEdge(i+1);
360 bx[nx] = binsx->GetXmax();
361 by[ny] = binsy->GetXmax();
362 bz[nz] = binsz->GetXmax();
363
364 SetBinning(h, bx, by, bz);
365}
366
367// --------------------------------------------------------------------------
368//
369// Applies the binnings of one root-histogram x to another one h
370// Both histograms must be of the same type: TH1, TH2 or TH3
371//
372void MH::SetBinning(TH1 *h, const TH1 *x)
373{
374 if (h->InheritsFrom(TH3::Class()) && x->InheritsFrom(TH3::Class()))
375 {
376 SetBinning((TH3*)h, ((TH1*)x)->GetXaxis(), ((TH1*)x)->GetYaxis(), ((TH1*)x)->GetZaxis());
377 return;
378 }
379 if (h->InheritsFrom(TH3::Class()) || x->InheritsFrom(TH3::Class()))
380 return;
381 if (h->InheritsFrom(TH2::Class()) && x->InheritsFrom(TH2::Class()))
382 {
383 SetBinning((TH2*)h, ((TH1*)x)->GetXaxis(), ((TH1*)x)->GetYaxis());
384 return;
385 }
386 if (h->InheritsFrom(TH2::Class()) || x->InheritsFrom(TH2::Class()))
387 return;
388 if (h->InheritsFrom(TH1::Class()) && x->InheritsFrom(TH1::Class()))
389 {
390 SetBinning(h, ((TH1*)x)->GetXaxis());
391 return;
392 }
393}
394
395// --------------------------------------------------------------------------
396//
397// Multiplies all entries in a TArrayD by a float f
398//
399void MH::ScaleArray(TArrayD &bins, Double_t f)
400{
401 for (int i=0; i<bins.GetSize(); i++)
402 bins[i] *= f;
403}
404
405// --------------------------------------------------------------------------
406//
407// Scales the binning of a TAxis by a float f
408//
409TArrayD MH::ScaleAxis(TAxis &axe, Double_t f)
410{
411 TArrayD arr(axe.GetNbins()+1);
412
413 for (int i=1; i<=axe.GetNbins()+1; i++)
414 arr[i-1] = axe.GetBinLowEdge(i);
415
416 ScaleArray(arr, f);
417
418 return arr;
419}
420
421// --------------------------------------------------------------------------
422//
423// Scales the binning of one, two or three axis of a histogram by a float f
424//
425void MH::ScaleAxis(TH1 *h, Double_t fx, Double_t fy, Double_t fz)
426{
427 if (h->InheritsFrom(TH3::Class()))
428 {
429 SetBinning((TH3*)h,
430 ScaleAxis(*h->GetXaxis(), fx),
431 ScaleAxis(*h->GetYaxis(), fy),
432 ScaleAxis(*h->GetZaxis(), fz));
433 return;
434 }
435
436 if (h->InheritsFrom(TH2::Class()))
437 {
438 SetBinning((TH2*)h,
439 ScaleAxis(*h->GetXaxis(), fx),
440 ScaleAxis(*h->GetYaxis(), fy));
441 return;
442 }
443
444 if (h->InheritsFrom(TH1::Class()))
445 SetBinning(h, ScaleAxis(*h->GetXaxis(), fx));
446}
447
448// --------------------------------------------------------------------------
449//
450// Tries to find a MBinning container with the name "Binning"+name
451// in the given parameter list. If it was found it is applied to the
452// given histogram. This is only valid for 1D-histograms
453//
454Bool_t MH::ApplyBinning(const MParList &plist, TString name, TH1 *h)
455{
456 if (h->InheritsFrom(TH2::Class()) || h->InheritsFrom(TH3::Class()))
457 {
458 gLog << warn << "MH::ApplyBinning: '" << h->GetName() << "' is not a basic TH1 object... no binning applied." << endl;
459 return kFALSE;
460 }
461
462 const MBinning *bins = (MBinning*)plist.FindObject("Binning"+name);
463 if (!bins)
464 {
465 gLog << warn << "Object 'Binning" << name << "' [MBinning] not found... no binning applied." << endl;
466 return kFALSE;
467 }
468
469 SetBinning(h, bins);
470 return kTRUE;
471}
472
473void MH::FindGoodLimits(Int_t nbins, Int_t &newbins, Double_t &xmin, Double_t &xmax, Bool_t isInteger)
474{
475#if ROOT_VERSION_CODE > ROOT_VERSION(3,04,01)
476 THLimitsFinder::OptimizeLimits(nbins, newbins, xmin, xmax, isInteger);
477#else
478//*-*-*-*-*-*-*-*-*Find reasonable bin values*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
479//*-* ==========================
480
481 Double_t dx = 0.1*(xmax-xmin);
482 Double_t umin = xmin - dx;
483 Double_t umax = xmax + dx;
484
485 if (umin < 0 && xmin >= 0)
486 umin = 0;
487
488 if (umax > 0 && xmax <= 0)
489 umax = 0;
490
491 Double_t binlow =0;
492 Double_t binhigh =0;
493 Double_t binwidth=0;
494
495 TGaxis::Optimize(umin, umax, nbins, binlow, binhigh, nbins, binwidth, "");
496
497 if (binwidth <= 0 || binwidth > 1.e+39)
498 {
499 xmin = -1;
500 xmax = 1;
501 }
502 else
503 {
504 xmin = binlow;
505 xmax = binhigh;
506 }
507
508 if (isInteger)
509 {
510 Int_t ixmin = (Int_t)xmin;
511 Int_t ixmax = (Int_t)xmax;
512 Double_t dxmin = (Double_t)ixmin;
513 Double_t dxmax = (Double_t)ixmax;
514
515 xmin = xmin<0 && xmin!=dxmin ? dxmin - 1 : dxmin;
516 xmax = xmax>0 && xmax!=dxmax ? dxmax + 1 : dxmax;
517
518 if (xmin>=xmax)
519 xmax = xmin+1;
520
521 Int_t bw = 1 + (Int_t)((xmax-xmin)/nbins);
522
523 nbins = (Int_t)((xmax-xmin)/bw);
524
525 if (xmin+nbins*bw < xmax)
526 {
527 nbins++;
528 xmax = xmin +nbins*bw;
529 }
530 }
531
532 newbins = nbins;
533#endif
534}
535
536// --------------------------------------------------------------------------
537//
538// Returns the lowest entry in a histogram which is greater than gt (eg >0)
539//
540Double_t MH::GetMinimumGT(const TH1 &h, Double_t gt)
541{
542 Double_t min = FLT_MAX;
543
544 const TAxis &axe = *((TH1&)h).GetXaxis();
545
546 for (int i=1; i<=axe.GetNbins(); i++)
547 {
548 Double_t x = h.GetBinContent(i);
549 if (gt<x && x<min)
550 min = x;
551 }
552 return min;
553}
554
555// --------------------------------------------------------------------------
556//
557// Returns the bin center in a logarithmic scale. If the given bin
558// number is <1 it is set to 1. If it is =GetNbins() it is set to
559// GetNbins()
560//
561Double_t MH::GetBinCenterLog(const TAxis &axe, Int_t nbin)
562{
563 if (nbin>axe.GetNbins())
564 nbin = axe.GetNbins();
565
566 if (nbin<1)
567 nbin = 1;
568
569 const Double_t lo = axe.GetBinLowEdge(nbin);
570 const Double_t hi = axe.GetBinUpEdge(nbin);
571
572 const Double_t val = log10(lo) + log10(hi);
573
574 return pow(10, val/2);
575}
576
577// --------------------------------------------------------------------------
578//
579// Draws a copy of the two given histograms. Uses title as the pad title.
580// Also layout the two statistic boxes and a legend.
581//
582void MH::DrawCopy(const TH1 &hist1, const TH1 &hist2, const TString title)
583{
584 //
585 // Draw first histogram
586 //
587 TH1 *h1 = ((TH1&)hist1).DrawCopy();
588 gPad->SetBorderMode(0);
589 gPad->Update();
590
591 // FIXME: Also align max/min with set Maximum/Minimum
592 const Double_t maxbin1 = hist1.GetBinContent(hist1.GetMaximumBin());
593 const Double_t maxbin2 = hist2.GetBinContent(hist2.GetMaximumBin());
594 const Double_t minbin1 = hist1.GetBinContent(hist1.GetMinimumBin());
595 const Double_t minbin2 = hist2.GetBinContent(hist2.GetMinimumBin());
596
597 const Double_t max = TMath::Max(maxbin1, maxbin2);
598 const Double_t min = TMath::Min(minbin1, minbin2);
599
600 h1->SetMaximum(max>0?max*1.05:max*0.95);
601 h1->SetMinimum(max>0?min*0.95:min*1.05);
602
603 TPaveText *t = (TPaveText*)gPad->FindObject("title");
604 if (t)
605 {
606 t->SetName((TString)"MHTitle"); // rename object
607 t->Clear(); // clear old lines
608 t->AddText((TString)" "+title+" "); // add the new title
609 t->SetBit(kCanDelete); // make sure object is deleted
610
611 //
612 // FIXME: This is a stupid workaround to hide the redrawn
613 // (see THistPainter::PaintTitle) title
614 //
615 gPad->Modified(); // indicates a change
616 gPad->Update(); // recreates the original title
617 t->Pop(); // bring our title on top
618 }
619
620 //
621 // Rename first statistics box
622 //
623 TPaveStats &s1 = *(TPaveStats*)gPad->FindObject("stats");
624 const Double_t x1 = s1.GetX1NDC()-0.01;
625 s1.SetName((TString)"Stat"+hist1.GetTitle());
626 s1.SetX1NDC(x1-(s1.GetX2NDC()-s1.GetX1NDC()));
627 s1.SetX2NDC(x1);
628
629 //
630 // Draw second histogram
631 //
632 ((TH1&)hist2).DrawCopy("sames");
633 gPad->Update();
634
635 //
636 // Draw Legend
637 //
638 TPaveStats &s2 = *(TPaveStats*)gPad->FindObject("stats");
639 TLegend &l = *new TLegend(s2.GetX1NDC(),
640 s2.GetY1NDC()-0.015-(s2.GetY2NDC()-s2.GetY1NDC())/2,
641 s2.GetX2NDC(),
642 s2.GetY1NDC()-0.01
643 );
644 l.AddEntry((TH1*)&hist1, hist1.GetTitle());
645 l.AddEntry((TH1*)&hist2, hist2.GetTitle());
646 l.SetTextSize(s2.GetTextSize());
647 l.SetTextFont(s2.GetTextFont());
648 l.SetBorderSize(s2.GetBorderSize());
649 l.SetBit(kCanDelete);
650 l.Draw();
651}
652
653// --------------------------------------------------------------------------
654//
655// Draws the two given histograms. Uses title as the pad title.
656// Also layout the two statistic boxes and a legend.
657//
658void MH::Draw(TH1 &hist1, TH1 &hist2, const TString title)
659{
660 //
661 // Draw first histogram
662 //
663 hist1.Draw();
664 gPad->SetBorderMode(0);
665 gPad->Update();
666
667 if (hist1.GetEntries()>0 && hist2.GetEntries()>0)
668 {
669 const Double_t maxbin1 = hist1.GetBinContent(hist1.GetMaximumBin());
670 const Double_t maxbin2 = hist2.GetBinContent(hist2.GetMaximumBin());
671 const Double_t minbin1 = hist1.GetBinContent(hist1.GetMinimumBin());
672 const Double_t minbin2 = hist2.GetBinContent(hist2.GetMinimumBin());
673
674 const Double_t max = TMath::Max(maxbin1, maxbin2);
675 const Double_t min = TMath::Min(minbin1, minbin2);
676
677 if (max!=min)
678 {
679 hist1.SetMaximum(max>0?max*1.05:max*0.95);
680 hist1.SetMinimum(max>0?min*0.95:min*1.05);
681 }
682 }
683
684 TPaveText *t = (TPaveText*)gPad->FindObject("title");
685 if (t)
686 {
687 t->SetName((TString)"MHTitle"); // rename object
688 t->Clear(); // clear old lines
689 t->AddText((TString)" "+title+" "); // add the new title
690 t->SetBit(kCanDelete); // make sure object is deleted
691
692 //
693 // FIXME: This is a stupid workaround to hide the redrawn
694 // (see THistPainter::PaintTitle) title
695 //
696 gPad->Modified(); // indicates a change
697 gPad->Update(); // recreates the original title
698 t->Pop(); // bring our title on top
699 }
700
701 //
702 // Rename first statistics box
703 //
704 TPaveStats &s1 = *(TPaveStats*)gPad->FindObject("stats");
705 const Double_t x1 = s1.GetX1NDC()-0.01;
706 s1.SetName((TString)"Stat"+hist1.GetTitle());
707 s1.SetX1NDC(x1-(s1.GetX2NDC()-s1.GetX1NDC()));
708 s1.SetX2NDC(x1);
709
710 //
711 // Draw second histogram
712 //
713 hist2.Draw("sames");
714 gPad->Update();
715
716 //
717 // Draw Legend
718 //
719 TPaveStats &s2 = *(TPaveStats*)gPad->FindObject("stats");
720 TLegend &l = *new TLegend(s2.GetX1NDC(),
721 s2.GetY1NDC()-0.015-(s2.GetY2NDC()-s2.GetY1NDC())/2,
722 s2.GetX2NDC(),
723 s2.GetY1NDC()-0.01
724 );
725 l.AddEntry(&hist1, hist1.GetTitle());
726 l.AddEntry(&hist2, hist2.GetTitle());
727 l.SetTextSize(s2.GetTextSize());
728 l.SetTextFont(s2.GetTextFont());
729 l.SetBorderSize(s2.GetBorderSize());
730 l.SetBit(kCanDelete);
731 l.Draw();
732}
733
734// --------------------------------------------------------------------------
735//
736// If the opt string contains 'nonew' or gPad is not given NULL is returned.
737// Otherwise the present gPad is returned.
738//
739TVirtualPad *MH::GetNewPad(Option_t *opt)
740{
741 TString str(opt);
742
743 if (!str.Contains("nonew", TString::kIgnoreCase))
744 return NULL;
745
746 return gPad;
747}
748
749// --------------------------------------------------------------------------
750//
751// If the opt string contains 'nonew' or gPad is not given a new canvas
752// with size w/h is created. Otherwise the object is cloned and drawn
753// to the present pad. The kCanDelete bit is set for the clone.
754//
755TObject *MH::DrawClone(Option_t *opt, Int_t w, Int_t h) const
756{
757 TVirtualPad *p = GetNewPad(opt);
758 if (!p)
759 p = MakeDefCanvas(this, w, h);
760 else
761 p->Clear();
762
763 gROOT->SetSelectedPad(NULL);
764
765 TObject *o = MParContainer::DrawClone(opt);
766 o->SetBit(kCanDelete);
767 return o;
768}
769
770// --------------------------------------------------------------------------
771//
772// Check whether a class inheriting from MH overwrites the Draw function
773//
774Bool_t MH::OverwritesDraw(TClass *cls) const
775{
776 if (!cls)
777 cls = IsA();
778
779 //
780 // Check whether we reached the base class MTask
781 //
782 if (TString(cls->GetName())=="MH")
783 return kFALSE;
784
785 //
786 // Check whether the class cls overwrites Draw
787 //
788 if (cls->GetMethodAny("Draw"))
789 {
790 *fLog << all << "FOUND: " << cls->GetName() << " " << (cls->GetName()=="MH") << endl;
791 return kTRUE;
792 }
793
794 //
795 // If the class itself doesn't overload it check all it's base classes
796 //
797 TBaseClass *base=NULL;
798 TIter NextBase(cls->GetListOfBases());
799 while ((base=(TBaseClass*)NextBase()))
800 {
801 if (OverwritesDraw(base->GetClassPointer()))
802 return kTRUE;
803 }
804
805 return kFALSE;
806}
Note: See TracBrowser for help on using the repository browser.