source: trunk/Mars/mhbase/MH3.cc@ 12509

Last change on this file since 12509 was 10283, checked in by tbretz, 14 years ago
Added asignment of conversion function of TProfile2D.
File size: 41.7 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 2002 <mailto:tbretz@astro.uni-wuerzburg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2010
21!
22!
23\* ======================================================================== */
24
25/////////////////////////////////////////////////////////////////////////////
26//
27// MH3
28//
29// With this histogram you can fill a histogram with up to three
30// variables from Mars parameter containers in an eventloop.
31//
32// In the constructor you can give up to three variables which should be
33// filled in the histogram. Dependend on the number of given variables
34// (data members) a TH1D, TH2D or TH3D is created.
35// Specify the data mamber like the following:
36// "MHillas.fLength"
37// Where MHillas is the name of the parameter container in the parameter
38// list and fLength is the name of the data member which should be filled
39// in the histogram. Assuming that your MHillas container has a different
40// name (MyHillas) the name to give would be:
41// "MyHillas.fLength"
42//
43// If you want to use a different unit for histogramming use SetScaleX,
44// SetScaleY and SetScaleZ.
45//
46//
47// Binning/Binning name
48// =====================
49//
50// The axis binning is retrieved from the parameter list, too. Create a
51// MBinning with the name "Binning" plus the name of your MH3 container
52// plus the axis name ("X", "Y" or "Z") and add it to the parameter list.
53//
54// If the binning should have a different name than the histogram name
55// the binning name can be added to the name, eg.:
56// SetName("MyHistName;MyXBins;MyYBins")
57// Instead of BinningMyHistName[XYZ] the parameter list will be searched
58// for BinningMyXBinning, BinningMyYBins and BinningMyHistNameZ
59//
60// If you don't want to use a MBinning object from the parameter list
61// you can also set one directly, for example
62// MBinning bins(10, 0, 1);
63// mh3.SetBinningX(&bins);
64// You must not delete the MBinning object before the class has been
65// PreProcessed.
66//
67//
68// Axis titles
69// ===========
70//
71// 1) If no other title is given the rule for this axis is used.
72// 2) If the MBinning used for this axis has a non-default Title
73// (MBinning::HasTitle) this title is used for the corresponding axis
74// 3) If the MH3 has a non-default title (MH3::SetTitle called)
75// this title is set as the histogram title. It can be used to overwrite
76// the axis titles. For more information see TH1::SetTitle, eg.
77// SetTitle("MyHist;x[mm];y[cm];Counts");
78//
79//
80// Labels
81// ======
82//
83// To use labels at an axis you have to initialize this for the axis
84// by either calling InitLabels(Labels_t) or setiting a DefaultLabel.
85// For the axis for which the labels have been initialized the
86// number returned by the given corresponding phrase is converted
87// to int with TMath::Nint and used as a label. If you want to replace
88// this id by a named label you can call DefineLabel to do that.
89// Several ids can be replaced by the same label. If you define
90// named labels for every label which was not defined the default
91// is used, if any, otherwise an unnamed label is created.
92//
93// In the case of an axis with labels the axis-title cannot be
94// set via a MBinning, because the MBinning is not evaluated.
95//
96// Please note that for some reason not all combinations of
97// labels, dimensions and weights are available in the root-
98// histogram classes. Please check the MH3::Fill function to see
99// whether your combination is supported.
100//
101//
102// Examples:
103// =========
104//
105// 1) MH3 myhist("MHillas.fLength");
106// myhist.SetName("MyHist");
107// myhist.SetScaleX(geomcam.GetConvMm2Deg()); //convert length to degree
108// MBinning bins("BinningMyHistX", "Title for my x-axis [Hz]");
109// bins.SetEdges(10, 0, 150);
110// plist.AddToList(&bins);
111//
112// 2) MH3 myhist("MHillas.fLength");
113// myhist.SetName("MyHist;MyX");
114// myhist.SetTitle("Histogram Title;X-Title [mm];Counts");
115// MBinning bins("BinningMyX");
116// bins.SetEdges(10, 0, 150);
117// plist.AddToList(&bins);
118//
119// 3) MH3 myhist("MTriggerPatter.GetUnprescaled");
120// myhist.SetWeight("1./MRawRunHeader.GetRunLength");
121// myhist.SetTitle("Rate of the trigger pattern [Hz];Run Number;Trigger Pattern;Rate [Hz]");
122// myhist.InitLabels(MH3::kLabelsXY);
123// myhist.DefaultLabelY("UNKNOWN"); // Lvl1
124// myhist.DefineLabelY( 1, "Trig"); // Lvl1
125// myhist.DefineLabelY( 2, "Cal"); // Cal
126// myhist.DefineLabelY( 4, "Trig"); // Lvl2
127// myhist.DefineLabelY( 8, "Ped"); // Ped
128//
129//
130// Class Version 1:
131// ----------------
132// - MData *fData[3];
133// + MData *fData[4];
134//
135// Class Version 2:
136// ----------------
137// - MDataChain *fData[3]; // Object from which the data is filled
138// + MData *fData[3]; // Object from which the data is filled
139//
140// Class Version 3:
141// ----------------
142// - Byte_t fStyleBits
143// + MBinning fBins[3]
144//
145// Class Version 5:
146// ----------------
147// + TFormula *fConversion
148//
149// Class Version 6:
150// ----------------
151// + MData *fWeight;
152//
153// Class Version 7:
154// ----------------
155// + MData *fData[4]
156// + Double_t fScale[4]
157// - MData *fData[3]
158// - Double_t fScale[3]
159//
160/////////////////////////////////////////////////////////////////////////////
161#include "MH3.h"
162
163#include <ctype.h> // tolower
164#include <stdlib.h> // atoi (Ubuntu 8.10)
165#include <fstream>
166
167#include <TMath.h>
168#include <TFormula.h>
169
170#include <THashList.h>
171#include <TObjString.h>
172
173//#include <TPad.h>
174#include <TStyle.h>
175#include <TCanvas.h>
176
177#include <TH2.h>
178#include <TH3.h>
179#include <TProfile.h>
180#include <TProfile2D.h>
181#include <TProfile3D.h>
182
183#include "MLog.h"
184#include "MLogManip.h"
185
186#include "MString.h"
187
188#include "MParList.h"
189#include "MBinning.h"
190#include "MDataPhrase.h"
191
192ClassImp(MH3);
193
194using namespace std;
195
196const TString MH3::gsDefName = "MH3";
197const TString MH3::gsDefTitle = "Container for a n-D Mars Histogram";
198
199// --------------------------------------------------------------------------
200//
201// Set fStyleBits to 0, Reset fBins to NULL, fScale to 1, set name and title
202// to gsDefName and gsDefTitle and if fHist!=NULL UseCurrentStyle and
203// SetDirectory(0)
204//
205void MH3::Init()
206{
207 fStyleBits = 0;
208
209 fWeight = NULL;
210
211 fData[0] = NULL;
212 fData[1] = NULL;
213 fData[2] = NULL;
214 fData[3] = NULL;
215
216 fBins[0] = NULL;
217 fBins[1] = NULL;
218 fBins[2] = NULL;
219
220 fScale[0] = 1;
221 fScale[1] = 1;
222 fScale[2] = 1;
223 fScale[3] = 1;
224
225 fConversion = NULL;
226
227 fName = gsDefName;
228 fTitle = gsDefTitle;
229
230 if (!fHist)
231 return;
232
233 fHist->UseCurrentStyle();
234 fHist->SetDirectory(NULL);
235}
236
237// --------------------------------------------------------------------------
238//
239// Default constructor.
240//
241MH3::MH3(const Int_t dim, Type_t type) : fDimension(dim), fHist(NULL)
242{
243 // FIXME?
244 switch (type)
245 {
246 case kHistogram:
247 if (fDimension>3)
248 fDimension=3;
249 break;
250 case kProfile:
251 case kProfileSpread:
252 fDimension = -TMath::Abs(fDimension);
253 if (fDimension<-2)
254 fDimension = -2;
255 break;
256 }
257
258 switch (fDimension)
259 {
260 case 1:
261 fHist = new TH1D;
262 fHist->SetYTitle("Counts");
263 break;
264 case -1:
265 fHist = new TProfile;
266 fHist->SetYTitle("Average");
267 if (type==kProfileSpread)
268 static_cast<TProfile*>(fHist)->SetErrorOption("s");
269 break;
270 case 2:
271 fHist = new TH2D;
272 fHist->SetZTitle("Counts");
273 break;
274 case -2:
275 fHist = new TProfile2D;
276 fHist->SetZTitle("Average");
277 static_cast<TProfile2D*>(fHist)->BuildOptions(0, 0, type==kProfileSpread?"s":"");
278 break;
279 case 3:
280 fHist = new TH3D;
281 break;
282 case -3:
283 fHist = new TProfile2D;
284 fHist->SetZTitle("Average");
285 static_cast<TProfile2D*>(fHist)->BuildOptions(0, 0, type==kProfileSpread?"s":"");
286 break;
287 }
288
289 Init();
290}
291
292// --------------------------------------------------------------------------
293//
294// Creates an TH1D. memberx is filled into the X-bins. For a more detailed
295// description see the class description above.
296//
297MH3::MH3(const char *memberx, Type_t type) : fDimension(1)
298{
299 fHist = new TH1D;
300 fHist->SetYTitle("Counts");
301
302 Init();
303
304 fData[0] = new MDataPhrase(memberx);
305}
306
307// --------------------------------------------------------------------------
308//
309// Adapt a given histogram
310//
311MH3::MH3(const TH1 &h1) : fDimension(1)
312{
313 if (h1.InheritsFrom(TH3::Class()))
314 fDimension = 3;
315 if (h1.InheritsFrom(TH2::Class()))
316 fDimension = 2;
317
318 if (h1.InheritsFrom(TProfile3D::Class()))
319 {
320 fDimension = -3;
321 *fLog << warn << "WARNING - MH3::MH3(TH1&) does not support TProfile3D." << endl;
322 }
323 if (h1.InheritsFrom(TProfile2D::Class()))
324 fDimension = -2;
325 if (h1.InheritsFrom(TProfile::Class()))
326 fDimension = -1;
327
328 fHist = (TH1*)h1.Clone();
329
330 Init(); // Before without SeUseCurrentStyle!
331
332 switch (fDimension)
333 {
334 case 3:
335 case -2:
336 fData[2] = new MDataPhrase(h1.GetZaxis()->GetTitle());
337 case 2:
338 case -1:
339 fData[1] = new MDataPhrase(h1.GetYaxis()->GetTitle());
340 case 1:
341 fData[0] = new MDataPhrase(h1.GetXaxis()->GetTitle());
342 }
343}
344
345// --------------------------------------------------------------------------
346//
347// Creates an TH2D. memberx is filled into the X-bins. membery is filled
348// into the Y-bins. For a more detailed description see the class
349// description above.
350//
351MH3::MH3(const char *memberx, const char *membery, Type_t type)
352 : fDimension(type==kHistogram?2:-1)
353{
354
355 switch (fDimension)
356 {
357 case 2:
358 fHist = static_cast<TH1*>(new TH2D);
359 break;
360 case -1:
361 fHist = static_cast<TH1*>(new TProfile);
362 if (type==kProfileSpread)
363 static_cast<TProfile*>(fHist)->SetErrorOption("s");
364
365 break;
366 }
367
368 fHist->SetZTitle(fDimension>0?"Counts":"Average");
369
370 Init();
371
372 fData[0] = new MDataPhrase(memberx);
373 fData[1] = new MDataPhrase(membery);
374}
375
376// --------------------------------------------------------------------------
377//
378// Creates an TH3D. memberx is filled into the X-bins. membery is filled
379// into the Y-bins. membery is filled into the Z-bins. For a more detailed
380// description see the class description above.
381//
382MH3::MH3(const char *memberx, const char *membery, const char *memberz, Type_t type)
383 : fDimension(type==kHistogram?3:-2)
384{
385 switch (fDimension)
386 {
387 case 3:
388 fHist = static_cast<TH1*>(new TH3D);
389 break;
390 case -2:
391 fHist = static_cast<TH1*>(new TProfile2D);
392 static_cast<TProfile2D*>(fHist)->BuildOptions(0, 0, type==kProfileSpread?"s":"");
393 }
394
395 Init();
396
397 fData[0] = new MDataPhrase(memberx);
398 fData[1] = new MDataPhrase(membery);
399 fData[2] = new MDataPhrase(memberz);
400}
401
402// --------------------------------------------------------------------------
403//
404// Creates an TH3D. memberx is filled into the X-bins. membery is filled
405// into the Y-bins. membery is filled into the Z-bins. Weight is used as a
406// weight for the profile histogram. For a more detailed description see the
407// class description above.
408//
409MH3::MH3(const char *memberx, const char *membery, const char *memberz, const char *weight, Type_t type)
410 : fDimension(-3)
411{
412 fHist = static_cast<TH1*>(new TProfile3D);
413 static_cast<TProfile3D*>(fHist)->BuildOptions(0, 0, type==kProfileSpread?"s":"");
414
415 Init();
416
417 fData[0] = new MDataPhrase(memberx);
418 fData[1] = new MDataPhrase(membery);
419 fData[2] = new MDataPhrase(memberz);
420 fData[3] = new MDataPhrase(weight);
421}
422
423// --------------------------------------------------------------------------
424//
425// Deletes the histogram
426//
427MH3::~MH3()
428{
429 if (fHist)
430 delete fHist;
431
432 if (fConversion)
433 delete fConversion;
434
435 if (fWeight)
436 delete fWeight;
437
438 for (int i=0; i<4; i++)
439 if (fData[i])
440 delete fData[i];
441
442 for (int i=0; i<3; i++)
443 if (fLabels[i].GetDefault())
444 delete fLabels[i].GetDefault();
445}
446
447// --------------------------------------------------------------------------
448//
449// You can set a weight as a phrase additionally to the one given
450// as an argument to Fill (most likely from MFillH). The two weights
451// are multiplied together.
452//
453void MH3::SetWeight(const char *phrase)
454{
455 if (fWeight)
456 delete fWeight;
457 fWeight = new MDataPhrase(phrase);
458}
459
460// --------------------------------------------------------------------------
461//
462// Set a function which is applied to the histogram before it is displayed.
463// Note, that it only effects the displayed histogram.
464//
465// e.g. SetConversion("sqrt(x)");
466//
467Bool_t MH3::SetConversion(const char *func)
468{
469 if (TString(func).IsNull())
470 {
471 delete fConversion;
472 fConversion = 0;
473 return kTRUE;
474 }
475
476 fConversion = new TFormula;
477
478 // Must have a name otherwise all axis labels disappear like a miracle
479 fConversion->SetName("ConversionFunction");
480 if (fConversion->Compile(func))
481 {
482 *fLog << err << dbginf << "Syntax Error: TFormula::Compile failed for " << func << endl;
483 delete fConversion;
484 fConversion = 0;
485 return kFALSE;
486 }
487
488 gROOT->GetListOfFunctions()->Remove(fConversion);
489
490 return kTRUE;
491}
492
493// --------------------------------------------------------------------------
494//
495// The axis label is centered and the labeling of the axis is initialized.
496//
497// This function must not be called after any label has been created!
498//
499void MH3::InitLabels(TAxis &x) const
500{
501 x.CenterTitle();
502 x.SetBinLabel(1, "");
503 x.LabelsOption("h"); // FIXME: Is "a" thread safe? (Paint and Fill?)
504 x.GetLabels()->Delete();
505}
506
507// --------------------------------------------------------------------------
508//
509// Depending on the bits set the InitLabels(TAxis&) function for
510// the corresponding axes are called. In any case the kCanRebin bit
511// is set.
512//
513// This function must not be called after any label has been created!
514//
515void MH3::InitLabels(Labels_t type) const
516{
517 if (!fHist)
518 return;
519
520 if (type&kLabelsX && fHist->GetXaxis())
521 InitLabels(*fHist->GetXaxis());
522
523 if (type&kLabelsY && fHist->GetYaxis())
524 InitLabels(*fHist->GetYaxis());
525
526 if (type&kLabelsZ && fHist->GetZaxis())
527 InitLabels(*fHist->GetZaxis());
528
529 if (type&kLabelsXYZ)
530 fHist->SetBit(TH1::kCanRebin);
531}
532
533// --------------------------------------------------------------------------
534//
535// Return the corresponding Labels_t describing for which axis
536// axis-labels are switched on.
537//
538MH3::Labels_t MH3::GetLabels() const
539{
540 UInt_t type = kNoLabels;
541 if (fHist->GetXaxis() && fHist->GetXaxis()->GetLabels())
542 type |= kLabelsX;
543 if (fHist->GetYaxis() && fHist->GetYaxis()->GetLabels())
544 type |= kLabelsY;
545 if (fHist->GetZaxis() && fHist->GetZaxis()->GetLabels())
546 type |= kLabelsZ;
547 return (Labels_t)type;
548}
549
550// --------------------------------------------------------------------------
551//
552// Calls the LabelsDeflate from the histogram for all axes.
553// LabelsDeflate will just do nothing if the axis has no labels
554// initialized.
555//
556void MH3::DeflateLabels() const
557{
558 fHist->LabelsDeflate("X");
559 fHist->LabelsDeflate("Y");
560 fHist->LabelsDeflate("Z");
561}
562
563// --------------------------------------------------------------------------
564//
565// Returns the named label corresponding to the given value
566// and the given axis. The names are defined with the
567// DefineLabel-functions. if no name is defined the value
568// is converted to a string with %d and TMath::Nint.
569// If names are defined, but not for the given value, the default
570// label is returned instead. If no default is defined the
571// %d-converted string is returned.
572//
573const char *MH3::GetLabel(Int_t axe, Double_t val) const
574{
575 const Int_t v = TMath::Nint(val);
576
577 if (fLabels[axe].GetSize())
578 {
579 const char *l = fLabels[axe].GetObjName(v);
580 if (l)
581 return l;
582 }
583
584 return Form("%d", v);
585}
586
587// --------------------------------------------------------------------------
588//
589// Return the data members used by the data chain to be used in
590// MTask::AddBranchToList
591//
592TString MH3::GetDataMember() const
593{
594 TString str=fData[0]->GetDataMember();
595
596 for (int i=1; i<4; i++)
597 if (fData[i])
598 {
599 str += ";";
600 str += fData[i]->GetDataMember();
601 }
602
603 return str;
604}
605
606// --------------------------------------------------------------------------
607//
608// Setup the Binning for the histograms automatically if the correct
609// instances of MBinning are found in the parameter list
610// For a more detailed description see class description above.
611//
612Bool_t MH3::SetupFill(const MParList *plist)
613{
614 // reset histogram (necessary if the same eventloop is run more than once)
615 if (!TestBit(kDoNotReset))
616 fHist->Reset();
617
618 // Tokenize name into name and binnings names
619 TObjArray *tok = fName.Tokenize(";");
620
621 const TString name = (*tok)[0] ? (*tok)[0]->GetName() : fName.Data();
622
623 TString bx = (*tok)[1] ? (*tok)[1]->GetName() : Form("%sX", name.Data());
624 TString by = (*tok)[2] ? (*tok)[2]->GetName() : Form("%sY", name.Data());
625 TString bz = (*tok)[3] ? (*tok)[3]->GetName() : Form("%sZ", name.Data());
626
627 bx.Prepend("Binning");
628 by.Prepend("Binning");
629 bz.Prepend("Binning");
630
631 delete tok;
632
633 MBinning *binsx = NULL;
634 MBinning *binsy = NULL;
635 MBinning *binsz = NULL;
636
637 const Labels_t labels = GetLabels();
638
639 switch (TMath::Abs(fDimension))
640 {
641 case 3:
642 if (fData[2])
643 fHist->SetZTitle(fData[2]->GetTitle());
644 if (!(labels&kLabelsZ))
645 {
646 binsz = fBins[2] ? fBins[2] : (MBinning*)plist->FindObject(bz, "MBinning");
647 if (!binsz)
648 {
649 *fLog << err << dbginf << "MBinning '" << bz << "' not found... aborting." << endl;
650 return kFALSE;
651 }
652 if (binsz->HasTitle())
653 fHist->SetZTitle(binsz->GetTitle());
654 if (binsz->IsLogarithmic())
655 fHist->SetBit(kIsLogz);
656 }
657 case 2:
658 if (fData[1])
659 fHist->SetYTitle(fData[1]->GetTitle());
660 if (!(labels&kLabelsY))
661 {
662 binsy = fBins[1] ? fBins[1] : (MBinning*)plist->FindObject(by, "MBinning");
663 if (!binsy)
664 {
665 *fLog << err << dbginf << "MBinning '" << by << "' not found... aborting." << endl;
666 return kFALSE;
667 }
668 if (binsy->HasTitle())
669 fHist->SetYTitle(binsy->GetTitle());
670 if (binsy->IsLogarithmic())
671 fHist->SetBit(kIsLogy);
672 }
673 case 1:
674 if (fData[0]!=NULL)
675 fHist->SetXTitle(fData[0]->GetTitle());
676 if (!(labels&kLabelsX))
677 {
678 binsx = fBins[0] ? fBins[0] : (MBinning*)plist->FindObject(bx, "MBinning");
679 if (!binsx)
680 {
681 if (fDimension==1)
682 binsx = (MBinning*)plist->FindObject("Binning"+fName, "MBinning");
683
684 if (!binsx)
685 {
686 *fLog << err << dbginf << "Neither '" << bx << "' nor 'Binning" << fName << "' found... aborting." << endl;
687 return kFALSE;
688 }
689 }
690 if (binsx->HasTitle())
691 fHist->SetXTitle(binsx->GetTitle());
692 if (binsx->IsLogarithmic())
693 fHist->SetBit(kIsLogx);
694 }
695 }
696
697 // PreProcess existing fData members
698 for (int i=0; i<4; i++)
699 if (fData[i] && !fData[i]->PreProcess(plist))
700 return kFALSE;
701
702 if (fWeight && !fWeight->PreProcess(plist))
703 return kFALSE;
704
705 TString title(fDimension>0?"Histogram":"Profile");
706 title += " for ";
707 title += name;
708 title += Form(" (%dD)", TMath::Abs(fDimension));
709
710 fHist->SetName(name);
711 fHist->SetTitle(fTitle==gsDefTitle ? title : fTitle);
712 fHist->SetDirectory(0);
713
714 // This is for the case we have set lables
715 const MBinning def(1, 0, 1);
716 if (!binsx)
717 binsx = const_cast<MBinning*>(&def);
718 if (!binsy)
719 binsy = const_cast<MBinning*>(&def);
720 if (!binsz)
721 binsz = const_cast<MBinning*>(&def);
722
723 // set binning
724 switch (TMath::Abs(fDimension))
725 {
726 case 1:
727 SetBinning(*fHist, *binsx);
728 return kTRUE;
729 case 2:
730 SetBinning(static_cast<TH2&>(*fHist), *binsx, *binsy);
731 return kTRUE;
732 case 3:
733 SetBinning(static_cast<TH3&>(*fHist), *binsx, *binsy, *binsz);
734 return kTRUE;
735 }
736
737 *fLog << err << "ERROR - MH3 has " << TMath::Abs(fDimension) << " dimensions!" << endl;
738 return kFALSE;
739}
740
741// --------------------------------------------------------------------------
742//
743// Set the name of the histogram ant the MH3 container
744//
745void MH3::SetName(const char *name)
746{
747 if (fHist)
748 {
749 if (gPad)
750 {
751 const TString pfx(MString::Format("%sProfX", fHist->GetName()));
752 const TString pfy(MString::Format("%sProfY", fHist->GetName()));
753
754 TProfile *p = 0;
755 if ((p=dynamic_cast<TProfile*>(gPad->FindObject(pfx))))
756 p->SetName(MString::Format("%sProfX", name));
757 if ((p=dynamic_cast<TProfile*>(gPad->FindObject(pfy))))
758 p->SetName(MString::Format("%sProfY", name));
759 }
760
761 fHist->SetName(name);
762 fHist->SetDirectory(0);
763
764 }
765 MParContainer::SetName(name);
766}
767
768// --------------------------------------------------------------------------
769//
770// Set the title of the histogram ant the MH3 container
771//
772void MH3::SetTitle(const char *title)
773{
774 if (fHist)
775 fHist->SetTitle(title);
776 MParContainer::SetTitle(title);
777}
778
779// --------------------------------------------------------------------------
780//
781// Fills the one, two or three data members into our histogram
782//
783Int_t MH3::Fill(const MParContainer *par, const Stat_t ww)
784{
785 // Get Information about labels (UInt_t, to supress warning about
786 // unhandeled cases in switch)
787 const UInt_t type = GetLabels();
788
789 // Get values for axis
790 Double_t x=0;
791 Double_t y=0;
792 Double_t z=0;
793 Double_t t=0;
794 Double_t w=ww;
795
796 switch (fDimension)
797 {
798 case -3:
799 t = fData[3]->GetValue()*fScale[3];
800 case -2:
801 case 3:
802 z = fData[2]->GetValue()*fScale[2];
803 case -1:
804 case 2:
805 y = fData[1]->GetValue()*fScale[1];
806 case 1:
807 x = fData[0]->GetValue()*fScale[0];
808 }
809
810 if (fWeight)
811 w *= fWeight->GetValue();
812
813 // If label option is set, convert value to label
814 TString labelx, labely, labelz;
815 if (type&kLabelsX)
816 labelx = GetLabel(0, x);
817 if (type&kLabelsY)
818 labely = GetLabel(1, y);
819 if (type&kLabelsZ)
820 labelz = GetLabel(2, z);
821
822 // Fill histogram
823 switch (fDimension)
824 {
825 case 3:
826 switch (type)
827 {
828 case kNoLabels:
829 static_cast<TH3*>(fHist)->Fill(x, y, z, w);
830 return kTRUE;
831 case kLabelsX:
832 static_cast<TH3*>(fHist)->Fill(labelx, y, z);
833 return kTRUE;
834 case kLabelsY:
835 static_cast<TH3*>(fHist)->Fill(x, labely, z, w);
836 return kTRUE;
837 case kLabelsZ:
838 static_cast<TH3*>(fHist)->Fill(x, y, labelz, w);
839 return kTRUE;
840 case kLabelsXY:
841 static_cast<TH3*>(fHist)->Fill(labelx, labely, z, w);
842 return kTRUE;
843 case kLabelsXZ:
844 static_cast<TH3*>(fHist)->Fill(labelx, y, labelz, w);
845 return kTRUE;
846 case kLabelsYZ:
847 static_cast<TH3*>(fHist)->Fill(x, labely, labelz, w);
848 return kTRUE;
849 case kLabelsXYZ:
850 static_cast<TH3*>(fHist)->Fill(labelx, labely, labelz, w);
851 return kTRUE;
852 }
853 break;
854 case 2:
855 switch (type)
856 {
857 case kNoLabels:
858 static_cast<TH2*>(fHist)->Fill(x, y, w);
859 return kTRUE;
860 case kLabelsX:
861 static_cast<TH2*>(fHist)->Fill(x, labely, w);
862 return kTRUE;
863 case kLabelsY:
864 static_cast<TH2*>(fHist)->Fill(labelx, y, w);
865 return kTRUE;
866 case kLabelsXY:
867 static_cast<TH2*>(fHist)->Fill(labelx, labely, w);
868 return kTRUE;
869 }
870 break;
871 case 1:
872 switch (type)
873 {
874 case kNoLabels:
875 fHist->Fill(x, w);
876 return kTRUE;
877 case kLabelsX:
878 fHist->Fill(labelx, w);
879 return kTRUE;
880 }
881 break;
882 case -1:
883 switch (type)
884 {
885 case kNoLabels:
886 static_cast<TProfile*>(fHist)->Fill(x, y, w);
887 return kTRUE;
888 case kLabelsX:
889 static_cast<TProfile*>(fHist)->Fill(labelx, y, w);
890 return kTRUE;
891 }
892 break;
893 case -2:
894 switch (type)
895 {
896 case kNoLabels:
897 static_cast<TProfile2D*>(fHist)->Fill(x, y, z, w);
898 return kTRUE;
899 case kLabelsX:
900 static_cast<TProfile2D*>(fHist)->Fill(labelx, y, z);
901 return kTRUE;
902 case kLabelsY:
903 static_cast<TProfile2D*>(fHist)->Fill(x, labely, z);
904 return kTRUE;
905 case kLabelsXY:
906 static_cast<TProfile2D*>(fHist)->Fill(labelx, labely, z);
907 return kTRUE;
908 }
909 break;
910 case -3:
911 switch (type)
912 {
913 case kNoLabels:
914 static_cast<TProfile3D*>(fHist)->Fill(x, y, z, t, w);
915 return kTRUE;
916 default:
917 *fLog << err << "ERROR - Labels not supported in TProfile3D." << endl;
918 }
919 break;
920 }
921
922 *fLog << err << "MH3::Fill: ERROR - A fatal error occured." << endl;
923 return kERROR;
924}
925
926// --------------------------------------------------------------------------
927//
928// If an auto range bit is set the histogram range of the corresponding
929// axis is set to show only the non-empty bins (with a single empty bin
930// on both sides)
931//
932Bool_t MH3::Finalize()
933{
934 DeflateLabels();
935
936 Bool_t autorangex=TESTBIT(fStyleBits, 0);
937 Bool_t autorangey=TESTBIT(fStyleBits, 1);
938 //Bool_t autorangez=TESTBIT(fStyleBits, 2);
939
940 Int_t lo, hi;
941 if (autorangex)
942 {
943 GetRangeX(*fHist, lo, hi);
944 fHist->GetXaxis()->SetRange(lo-2, hi+1);
945 }
946 if (autorangey)
947 {
948 GetRangeY(*fHist, lo, hi);
949 fHist->GetYaxis()->SetRange(lo-2, hi+1);
950 }
951 /*
952 if (autorangez)
953 {
954 GetRangeZ(*fHist, lo, hi);
955 fHist->GetZaxis()->SetRange(lo-2, hi+1);
956 }
957 */
958 return kTRUE;
959}
960
961// --------------------------------------------------------------------------
962//
963// Apply the conversion function to the contents stored in fHist and
964// store the result in h.
965//
966// In the case of a TProfile we keep it a TProfile to keep the mean and
967// rms in y displayed. To get the error correctly we have to reverse
968// the calculation done in TProfile::GetBinError of course.
969//
970void MH3::Convert(TH1 &h) const
971{
972 const Bool_t prof = h.InheritsFrom(TProfile::Class()) || h.InheritsFrom(TProfile2D::Class()) || h.InheritsFrom(TProfile3D::Class());
973
974 for (Int_t z=0; z<=h.GetNbinsZ()+1; z++)
975 for (Int_t y=0; y<=h.GetNbinsY()+1; y++)
976 for (Int_t x=0; x<=h.GetNbinsX()+1; x++)
977 {
978 h.SetBinContent(x, y, z, fConversion->Eval(fHist->GetBinContent(x, y, z)));
979
980 if (prof)
981 h.SetBinError(x, y, z, TMath::Hypot(fConversion->Eval(fHist->GetBinContent(x, y, z)),
982 fConversion->Eval(fHist->GetBinError( x, y, z))));
983 else
984 h.SetBinError(x, y, z, fConversion->Eval(fHist->GetBinError(x, y, z)));
985 }
986
987 TProfile *p1 = dynamic_cast<TProfile*>(fHist);
988 if (p1)
989 for (Int_t i=0; i<p1->GetSize(); i++)
990 static_cast<TProfile&>(h).SetBinEntries(i, p1->GetBinEntries(i)>0 ? 1 : 0);
991
992 TProfile2D *p2 = dynamic_cast<TProfile2D*>(fHist);
993 if (p2)
994 for (Int_t i=0; i<p2->GetSize(); i++)
995 static_cast<TProfile2D&>(h).SetBinEntries(i, p2->GetBinEntries(i)>0 ? 1 : 0);
996
997 TProfile3D *p3 = dynamic_cast<TProfile3D*>(fHist);
998 if (p3)
999 for (Int_t i=0; i<p3->GetSize(); i++)
1000 static_cast<TProfile3D&>(h).SetBinEntries(i, p3->GetBinEntries(i)>0 ? 1 : 0);
1001}
1002
1003// --------------------------------------------------------------------------
1004//
1005// FIXME
1006//
1007void MH3::Paint(Option_t *o)
1008{
1009 TProfile *p=0;
1010
1011 if (TMath::Abs(fDimension)==2)
1012 MH::SetPalette("pretty");
1013
1014 if (fConversion)
1015 {
1016 TH1 *h = 0;
1017 if ((h=dynamic_cast<TH1*>(gPad->FindObject(fHist->GetName()))))
1018 Convert(*h);
1019 }
1020
1021 const TString pfx(MString::Format("%sProfX", fHist->GetName()));
1022 if ((p=dynamic_cast<TProfile*>(gPad->FindObject(pfx))))
1023 {
1024 Int_t col = p->GetLineColor();
1025 p = ((TH2*)fHist)->ProfileX(pfx, -1, -1, "s");
1026 p->SetLineColor(col);
1027 }
1028
1029 const TString pfy(MString::Format("%sProfY", fHist->GetName()));
1030 if ((p=dynamic_cast<TProfile*>(gPad->FindObject(pfy))))
1031 {
1032 Int_t col = p->GetLineColor();
1033 p = ((TH2*)fHist)->ProfileY(pfy, -1, -1, "s");
1034 p->SetLineColor(col);
1035 }
1036/*
1037 if (fHist->TestBit(kIsLogx) && fHist->GetEntries()>0) gPad->SetLogx();
1038 if (fHist->TestBit(kIsLogy) && fHist->GetEntries()>0) gPad->SetLogy();
1039 if (fHist->TestBit(kIsLogz) && fHist->GetEntries()>0) gPad->SetLogz();
1040 */
1041}
1042
1043// --------------------------------------------------------------------------
1044//
1045// If Xmax is < 3000*Xmin SetMoreLogLabels is called. If Xmax<5000
1046// the exponent is switched off (SetNoExponent)
1047//
1048void MH3::HandleLogAxis(TAxis &axe) const
1049{
1050 if (axe.GetXmax()>3000*axe.GetXmin())
1051 return;
1052
1053 axe.SetMoreLogLabels();
1054 if (axe.GetXmax()<5000)
1055 axe.SetNoExponent();
1056}
1057
1058// --------------------------------------------------------------------------
1059//
1060// Creates a new canvas and draws the histogram into it.
1061//
1062// Possible options are:
1063// PROFX: Draw a x-profile into the histogram (for 2D histograms only)
1064// PROFY: Draw a y-profile into the histogram (for 2D histograms only)
1065// ONLY: Draw the profile histogram only (for 2D histograms only)
1066// BLUE: Draw the profile in blue color instead of the histograms
1067// line color
1068//
1069// If the kIsLog?-Bit is set the axis is displayed lkogarithmically.
1070// eg this is set when applying a logarithmic MBinning
1071//
1072// Be careful: The histogram belongs to this object and won't get deleted
1073// together with the canvas.
1074//
1075void MH3::Draw(Option_t *opt)
1076{
1077 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(fHist);
1078 pad->SetBorderMode(0);
1079 pad->SetGridx();
1080 pad->SetGridy();
1081
1082 if (fHist->TestBit(kIsLogx))
1083 {
1084 pad->SetLogx();
1085 HandleLogAxis(*fHist->GetXaxis());
1086 }
1087 if (fHist->TestBit(kIsLogy))
1088 {
1089 pad->SetLogy();
1090 HandleLogAxis(*fHist->GetYaxis());
1091 }
1092 if (fHist->TestBit(kIsLogz))
1093 {
1094 pad->SetLogz();
1095 HandleLogAxis(*fHist->GetZaxis());
1096 }
1097
1098 fHist->SetFillStyle(4000);
1099
1100 TString str(opt);
1101 str.ToLower();
1102
1103 if (GetLabels())
1104 {
1105 if (str.IsNull() && fDimension==2)
1106 str = "colz";
1107
1108 if (str.Contains("box", TString::kIgnoreCase) && fDimension==2)
1109 fHist->SetLineColor(kBlue);
1110 }
1111
1112 const Bool_t only = str.Contains("only") && TMath::Abs(fDimension)==2;
1113 const Bool_t same = str.Contains("same") && TMath::Abs(fDimension)<3;
1114 const Bool_t blue = str.Contains("blue") && TMath::Abs(fDimension)==2;
1115 const Bool_t profx = str.Contains("profx") && TMath::Abs(fDimension)==2;
1116 const Bool_t profy = str.Contains("profy") && TMath::Abs(fDimension)==2;
1117
1118 str.ReplaceAll("only", "");
1119 str.ReplaceAll("blue", "");
1120 str.ReplaceAll("profx", "");
1121 str.ReplaceAll("profy", "");
1122 str.ReplaceAll(" ", "");
1123
1124 if (same && TMath::Abs(fDimension)==1)
1125 {
1126 fHist->SetLineColor(kBlue);
1127 fHist->SetMarkerColor(kBlue);
1128 }
1129
1130
1131 TH1 *h = fHist;
1132
1133 if (fConversion)
1134 {
1135 h = static_cast<TH1*>(fHist->Clone());
1136
1137 h->SetDirectory(0);
1138 h->SetBit(kCanDelete);
1139
1140 Convert(*h);
1141 }
1142
1143 // FIXME: We may have to remove all our own options from str!
1144 if (!only)
1145 h->Draw(str);
1146
1147 AppendPad();
1148
1149 TProfile *p=0;
1150 if (profx)
1151 {
1152 const TString pfx(MString::Format("%sProfX", h->GetName()));
1153
1154 if (same && (p=dynamic_cast<TProfile*>(gPad->FindObject(pfx))))
1155 *fLog << warn << "TProfile " << pfx << " already in pad." << endl;
1156
1157 p = ((TH2*)h)->ProfileX(pfx, -1, -1, "s");
1158 p->UseCurrentStyle();
1159 p->SetLineColor(blue ? kBlue : h->GetLineColor());
1160 p->SetBit(kCanDelete);
1161 p->SetDirectory(NULL);
1162 p->SetXTitle(h->GetXaxis()->GetTitle());
1163 p->SetYTitle(h->GetYaxis()->GetTitle());
1164 p->Draw(only&&!same?"":"same");
1165 }
1166 if (profy)
1167 {
1168 const TString pfy(MString::Format("%sProfY", h->GetName()));
1169
1170 if (same && (p=dynamic_cast<TProfile*>(gPad->FindObject(pfy))))
1171 *fLog << warn << "TProfile " << pfy << " already in pad." << endl;
1172
1173 p = ((TH2*)h)->ProfileY(pfy, -1, -1, "s");
1174 p->UseCurrentStyle();
1175 p->SetLineColor(blue ? kBlue : h->GetLineColor());
1176 p->SetBit(kCanDelete);
1177 p->SetDirectory(NULL);
1178 p->SetYTitle(h->GetXaxis()->GetTitle());
1179 p->SetXTitle(h->GetYaxis()->GetTitle());
1180 p->Draw(only&&!same?"":"same");
1181 }
1182
1183 //AppendPad("log");
1184}
1185
1186// --------------------------------------------------------------------------
1187//
1188// Implementation of SavePrimitive. Used to write the call to a constructor
1189// to a macro. In the original root implementation it is used to write
1190// gui elements to a macro-file.
1191//
1192void MH3::StreamPrimitive(ostream &out) const
1193{
1194 TString name = GetUniqueName();
1195
1196 out << " MH3 " << name << "(\"";
1197 out << fData[0]->GetRule() << "\"";
1198 if (fDimension>1 || fDimension<0)
1199 out << ", \"" << fData[1]->GetRule() << "\"";
1200 if (fDimension>2 || fDimension<-1)
1201 out << ", \"" << fData[2]->GetRule() << "\"";
1202
1203 out << ");" << endl;
1204
1205 if (fName!=gsDefName)
1206 out << " " << name << ".SetName(\"" << fName << "\");" << endl;
1207
1208 if (fTitle!=gsDefTitle)
1209 out << " " << name << ".SetTitle(\"" << fTitle << "\");" << endl;
1210
1211 if (fWeight)
1212 out << " " << name << ".SetWeight(\"" << fWeight->GetRule() << "\");" << endl;
1213
1214 switch (fDimension)
1215 {
1216 case -3:
1217 if (fScale[3]!=1)
1218 out << " " << name << ".SetScaleT(" << fScale[3] << ");" << endl;
1219 case -2:
1220 case 3:
1221 if (fScale[2]!=1)
1222 out << " " << name << ".SetScaleZ(" << fScale[2] << ");" << endl;
1223 case -1:
1224 case 2:
1225 if (fScale[1]!=1)
1226 out << " " << name << ".SetScaleY(" << fScale[1] << ");" << endl;
1227 case 1:
1228 if (fScale[0]!=1)
1229 out << " " << name << ".SetScaleX(" << fScale[0] << ");" << endl;
1230 }
1231}
1232
1233MH3::Type_t MH3::GetType() const
1234{
1235 switch (fDimension)
1236 {
1237 case -1:
1238 return TString(static_cast<TProfile*>(fHist)->GetErrorOption())=="s" ? kProfileSpread : kProfile;
1239 case -2:
1240 return TString(static_cast<TProfile2D*>(fHist)->GetErrorOption())=="s" ? kProfileSpread : kProfile;
1241 case -3:
1242 return TString(static_cast<TProfile3D*>(fHist)->GetErrorOption())=="s" ? kProfileSpread : kProfile;
1243 }
1244 return kHistogram;
1245}
1246
1247// --------------------------------------------------------------------------
1248//
1249// Used to rebuild a MH3 object of the same type (data members,
1250// dimension, ...)
1251//
1252MParContainer *MH3::New() const
1253{
1254 // FIXME: TREAT THE NEW OPTIONS CORRECTLY (PROFILE, LABELS)
1255
1256 MH3 *h = NULL;
1257
1258 if (fData[0] == NULL)
1259 h=new MH3(fDimension);
1260 else
1261 switch (fDimension)
1262 {
1263 case 1:
1264 h=new MH3(fData[0]->GetRule());
1265 break;
1266 case 2:
1267 case -1:
1268 h=new MH3(fData[0]->GetRule(), fData[1]->GetRule(), GetType());
1269 break;
1270 case 3:
1271 case -2:
1272 h=new MH3(fData[0]->GetRule(), fData[1]->GetRule(), fData[2]->GetRule(), GetType());
1273 break;
1274 case -3:
1275 h=new MH3(fData[0]->GetRule(), fData[1]->GetRule(), fData[2]->GetRule(), fData[3]->GetRule());
1276 break;
1277 }
1278
1279 switch (fDimension)
1280 {
1281 case -3:
1282 h->SetScaleZ(fScale[3]);
1283 case -2:
1284 case 3:
1285 h->SetScaleZ(fScale[2]);
1286 case 2:
1287 case -1:
1288 h->SetScaleY(fScale[1]);
1289 case 1:
1290 h->SetScaleX(fScale[0]);
1291 }
1292
1293 if (fWeight)
1294 h->SetWeight(fWeight->GetRule());
1295
1296 return h;
1297}
1298
1299// --------------------------------------------------------------------------
1300//
1301// FIXME
1302//
1303TString MH3::GetRule(const Char_t axis) const
1304{
1305 switch (tolower(axis))
1306 {
1307 case 'x':
1308 case 'X':
1309 return fData[0] ? fData[0]->GetRule() : TString("");
1310 case 'y':
1311 case 'Y':
1312 return fData[1] ? fData[1]->GetRule() : TString("");
1313 case 'z':
1314 case 'Z':
1315 return fData[2] ? fData[2]->GetRule() : TString("");
1316 case 't':
1317 case 'T':
1318 return fData[3] ? fData[3]->GetRule() : TString("");
1319 case 'w':
1320 case 'W':
1321 return fWeight ? fWeight->GetRule() : TString("");
1322 default:
1323 return "<n/a>";
1324 }
1325}
1326
1327// --------------------------------------------------------------------------
1328//
1329// Returns the total number of bins in a histogram (excluding under- and
1330// overflow bins)
1331//
1332Int_t MH3::GetNbins() const
1333{
1334 Int_t num = 1;
1335
1336 switch (TMath::Abs(fDimension))
1337 {
1338 case 3:
1339 num *= fHist->GetNbinsZ()+2;
1340 case 2:
1341 num *= fHist->GetNbinsY()+2;
1342 case 1:
1343 num *= fHist->GetNbinsX()+2;
1344 }
1345
1346 return num;
1347}
1348
1349// --------------------------------------------------------------------------
1350//
1351// Returns the total number of bins in a histogram (excluding under- and
1352// overflow bins) Return -1 if bin is underflow or overflow
1353//
1354Int_t MH3::FindFixBin(Double_t x, Double_t y, Double_t z) const
1355{
1356 const TAxis &axex = *fHist->GetXaxis();
1357 const TAxis &axey = *fHist->GetYaxis();
1358 const TAxis &axez = *fHist->GetZaxis();
1359
1360 Int_t binz = 0;
1361 Int_t biny = 0;
1362 Int_t binx = 0;
1363
1364 switch (fDimension)
1365 {
1366 case 3:
1367 case -2:
1368 binz = axez.FindFixBin(z);
1369 if (binz>axez.GetLast() || binz<axez.GetFirst())
1370 return -1;
1371 case 2:
1372 case -1:
1373 biny = axey.FindFixBin(y);
1374 if (biny>axey.GetLast() || biny<axey.GetFirst())
1375 return -1;
1376 case 1:
1377 binx = axex.FindFixBin(x);
1378 if (binx<axex.GetFirst() || binx>axex.GetLast())
1379 return -1;
1380 }
1381
1382 const Int_t nx = fHist->GetNbinsX()+2;
1383 const Int_t ny = fHist->GetNbinsY()+2;
1384
1385 return binx + nx*(biny +ny*binz);
1386}
1387
1388// --------------------------------------------------------------------------
1389//
1390// Return the MObjLookup corresponding to the axis/character.
1391// Note that only lower-case charecters (x, y, z) are supported.
1392// If for the axis no labels were set, the corresponding
1393// InitLabels is called.
1394//
1395MObjLookup *MH3::GetLabels(char axe)
1396{
1397 if (!fHist)
1398 return 0;
1399
1400 TAxis *x = 0;
1401
1402 switch (axe)
1403 {
1404 case 'x':
1405 x = fHist->GetXaxis();
1406 break;
1407 case 'y':
1408 x = fHist->GetYaxis();
1409 break;
1410 case 'z':
1411 x = fHist->GetZaxis();
1412 break;
1413 }
1414
1415 if (!x)
1416 return 0;
1417
1418 const Int_t idx = axe-'x';
1419
1420 if (!x->GetLabels())
1421 switch (idx)
1422 {
1423 case 0:
1424 InitLabels(kLabelsX);
1425 break;
1426 case 1:
1427 InitLabels(kLabelsY);
1428 break;
1429 case 2:
1430 InitLabels(kLabelsZ);
1431 break;
1432 }
1433
1434 return &fLabels[idx];
1435}
1436
1437// --------------------------------------------------------------------------
1438//
1439// Set a default label which is used if no other is found in the list
1440// of labels. if a default was set already it is overwritten. If the
1441// axis has not yet been initialized to use labels it it now.
1442//
1443void MH3::DefaultLabel(char axe, const char *name)
1444{
1445 MObjLookup *arr = GetLabels(axe);
1446 if (!arr)
1447 return;
1448
1449 if (arr->GetDefault())
1450 {
1451 delete arr->GetDefault();
1452 arr->SetDefault(0);
1453 }
1454
1455 if (name)
1456 arr->SetDefault(new TObjString(name));
1457}
1458
1459// --------------------------------------------------------------------------
1460//
1461// Define a name for a label. More than one label can have the same
1462// name. If the axis has not yet been initialized to use labels
1463// it it now.
1464//
1465void MH3::DefineLabel(char axe, Int_t label, const char *name)
1466{
1467 MObjLookup *arr = GetLabels(axe);
1468
1469 if (!arr || !name)
1470 return;
1471
1472 if (arr->GetObj(label)!=arr->GetDefault())
1473 return;
1474
1475 arr->Add(label, name);
1476}
1477
1478// --------------------------------------------------------------------------
1479//
1480// Define names for labels, like
1481// 1=Trig;2=Cal;4=Ped;8=Lvl2
1482// More than one label can have the same name. If the axis has not
1483// yet been initialized to use labels it it now.
1484//
1485// A default cannot be set here. Use DefaultLabel instead.
1486//
1487void MH3::DefineLabels(char axe, const TString &labels)
1488{
1489 TObjArray *arr = labels.Tokenize(';');
1490
1491 for (int i=0; i<arr->GetEntries(); i++)
1492 {
1493 const char *s = (*arr)[0]->GetName();
1494 const char *v = strchr(s, '=');
1495
1496 if (v)
1497 DefineLabel(axe, atoi(s), v+1);
1498 }
1499
1500 delete arr;
1501}
1502
1503void MH3::RecursiveRemove(TObject *obj)
1504{
1505 if (obj==fHist)
1506 fHist = 0;
1507
1508 MH::RecursiveRemove(obj);
1509}
Note: See TracBrowser for help on using the repository browser.