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

Last change on this file since 20103 was 20096, checked in by tbretz, 4 years ago
root 6.24 requires some more includes, clasng 10 is a bit more picky about overload warinings and casting, gErrorMutex i deprecated -> so we can only lock our error handler again ourselves... well... mess ahead... thank you.
File size: 42.1 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Thomas Bretz 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 {
431 delete fHist;
432 fHist = 0;
433 }
434
435 if (fConversion)
436 {
437 delete fConversion;
438 fConversion = 0;
439 }
440
441 if (fWeight)
442 {
443 delete fWeight;
444 fWeight = 0;
445 }
446
447 for (int i=0; i<4; i++)
448 if (fData[i])
449 {
450 delete fData[i];
451 fData[i] = 0;
452 }
453
454 for (int i=0; i<3; i++)
455 if (fLabels[i].GetDefault())
456 delete fLabels[i].GetDefault();
457}
458
459// --------------------------------------------------------------------------
460//
461// You can set a weight as a phrase additionally to the one given
462// as an argument to Fill (most likely from MFillH). The two weights
463// are multiplied together.
464//
465void MH3::SetWeight(const char *phrase)
466{
467 if (fWeight)
468 delete fWeight;
469 fWeight = new MDataPhrase(phrase);
470}
471
472// --------------------------------------------------------------------------
473//
474// Set a function which is applied to the histogram before it is displayed.
475// Note, that it only effects the displayed histogram.
476//
477// e.g. SetConversion("sqrt(x)");
478//
479Bool_t MH3::SetConversion(const char *func)
480{
481 if (TString(func).IsNull())
482 {
483 delete fConversion;
484 fConversion = 0;
485 return kTRUE;
486 }
487
488 fConversion = new TFormula;
489
490 // Must have a name otherwise all axis labels disappear like a miracle
491 fConversion->SetName("ConversionFunction");
492 if (fConversion->Compile(func))
493 {
494 *fLog << err << dbginf << "Syntax Error: TFormula::Compile failed for " << func << endl;
495 delete fConversion;
496 fConversion = 0;
497 return kFALSE;
498 }
499
500 gROOT->GetListOfFunctions()->Remove(fConversion);
501
502 return kTRUE;
503}
504
505// --------------------------------------------------------------------------
506//
507// The axis label is centered and the labeling of the axis is initialized.
508//
509// This function must not be called after any label has been created!
510//
511void MH3::InitLabels(TAxis &x) const
512{
513 x.CenterTitle();
514 x.SetBinLabel(1, "");
515 x.LabelsOption("h"); // FIXME: Is "a" thread safe? (Paint and Fill?)
516 x.GetLabels()->Delete();
517}
518
519// --------------------------------------------------------------------------
520//
521// Depending on the bits set the InitLabels(TAxis&) function for
522// the corresponding axes are called. In any case the kCanRebin bit
523// is set.
524//
525// This function must not be called after any label has been created!
526//
527void MH3::InitLabels(Labels_t type) const
528{
529 if (!fHist)
530 return;
531
532 if (type&kLabelsX && fHist->GetXaxis())
533 InitLabels(*fHist->GetXaxis());
534
535 if (type&kLabelsY && fHist->GetYaxis())
536 InitLabels(*fHist->GetYaxis());
537
538 if (type&kLabelsZ && fHist->GetZaxis())
539 InitLabels(*fHist->GetZaxis());
540
541#if ROOT_VERSION_CODE < ROOT_VERSION(6,00,00)
542 if (type&kLabelsXYZ)
543 fHist->SetBit(TH1::kCanRebin);
544#else
545 // https://root.cern.ch/content/main-histogram-changes-root-6
546 if (type&kLabelsXYZ)
547 fHist->SetCanExtend(TH1::kAllAxes);
548#endif
549}
550
551// --------------------------------------------------------------------------
552//
553// Return the corresponding Labels_t describing for which axis
554// axis-labels are switched on.
555//
556MH3::Labels_t MH3::GetLabels() const
557{
558 UInt_t type = kNoLabels;
559 if (fHist->GetXaxis() && fHist->GetXaxis()->GetLabels())
560 type |= kLabelsX;
561 if (fHist->GetYaxis() && fHist->GetYaxis()->GetLabels())
562 type |= kLabelsY;
563 if (fHist->GetZaxis() && fHist->GetZaxis()->GetLabels())
564 type |= kLabelsZ;
565 return (Labels_t)type;
566}
567
568// --------------------------------------------------------------------------
569//
570// Calls the LabelsDeflate from the histogram for all axes.
571// LabelsDeflate will just do nothing if the axis has no labels
572// initialized.
573//
574void MH3::DeflateLabels() const
575{
576 fHist->LabelsDeflate("X");
577 fHist->LabelsDeflate("Y");
578 fHist->LabelsDeflate("Z");
579}
580
581// --------------------------------------------------------------------------
582//
583// Returns the named label corresponding to the given value
584// and the given axis. The names are defined with the
585// DefineLabel-functions. if no name is defined the value
586// is converted to a string with %d and TMath::Nint.
587// If names are defined, but not for the given value, the default
588// label is returned instead. If no default is defined the
589// %d-converted string is returned.
590//
591const char *MH3::GetLabel(Int_t axe, Double_t val) const
592{
593 const Int_t v = TMath::Nint(val);
594
595 if (fLabels[axe].GetSize())
596 {
597 const char *l = fLabels[axe].GetObjName(v);
598 if (l)
599 return l;
600 }
601
602 return Form("%d", v);
603}
604
605// --------------------------------------------------------------------------
606//
607// Return the data members used by the data chain to be used in
608// MTask::AddBranchToList
609//
610TString MH3::GetDataMember() const
611{
612 TString str=fData[0]->GetDataMember();
613
614 for (int i=1; i<4; i++)
615 if (fData[i])
616 {
617 str += ";";
618 str += fData[i]->GetDataMember();
619 }
620
621 return str;
622}
623
624// --------------------------------------------------------------------------
625//
626// Setup the Binning for the histograms automatically if the correct
627// instances of MBinning are found in the parameter list
628// For a more detailed description see class description above.
629//
630Bool_t MH3::SetupFill(const MParList *plist)
631{
632 // reset histogram (necessary if the same eventloop is run more than once)
633 if (!TestBit(kDoNotReset))
634 fHist->Reset();
635
636 // Tokenize name into name and binnings names
637 TObjArray *tok = fName.Tokenize(";");
638
639 const TString name = (*tok)[0] ? (*tok)[0]->GetName() : fName.Data();
640
641 TString bx = (*tok)[1] ? (*tok)[1]->GetName() : Form("%sX", name.Data());
642 TString by = (*tok)[2] ? (*tok)[2]->GetName() : Form("%sY", name.Data());
643 TString bz = (*tok)[3] ? (*tok)[3]->GetName() : Form("%sZ", name.Data());
644
645 bx.Prepend("Binning");
646 by.Prepend("Binning");
647 bz.Prepend("Binning");
648
649 delete tok;
650
651 MBinning *binsx = NULL;
652 MBinning *binsy = NULL;
653 MBinning *binsz = NULL;
654
655 const Labels_t labels = GetLabels();
656
657 switch (TMath::Abs(fDimension))
658 {
659 case 3:
660 if (fData[2])
661 fHist->SetZTitle(fData[2]->GetTitle());
662 if (!(labels&kLabelsZ))
663 {
664 binsz = fBins[2] ? fBins[2] : (MBinning*)plist->FindObject(bz, "MBinning");
665 if (!binsz)
666 {
667 *fLog << err << dbginf << "MBinning '" << bz << "' not found... aborting." << endl;
668 return kFALSE;
669 }
670 if (binsz->HasTitle())
671 fHist->SetZTitle(binsz->GetTitle());
672 if (binsz->IsLogarithmic())
673 fHist->SetBit(kIsLogz);
674 }
675 case 2:
676 if (fData[1])
677 fHist->SetYTitle(fData[1]->GetTitle());
678 if (!(labels&kLabelsY))
679 {
680 binsy = fBins[1] ? fBins[1] : (MBinning*)plist->FindObject(by, "MBinning");
681 if (!binsy)
682 {
683 *fLog << err << dbginf << "MBinning '" << by << "' not found... aborting." << endl;
684 return kFALSE;
685 }
686 if (binsy->HasTitle())
687 fHist->SetYTitle(binsy->GetTitle());
688 if (binsy->IsLogarithmic())
689 fHist->SetBit(kIsLogy);
690 }
691 case 1:
692 if (fData[0]!=NULL)
693 fHist->SetXTitle(fData[0]->GetTitle());
694 if (!(labels&kLabelsX))
695 {
696 binsx = fBins[0] ? fBins[0] : (MBinning*)plist->FindObject(bx, "MBinning");
697 if (!binsx)
698 {
699 if (fDimension==1)
700 binsx = (MBinning*)plist->FindObject("Binning"+fName, "MBinning");
701
702 if (!binsx)
703 {
704 *fLog << err << dbginf << "Neither '" << bx << "' nor 'Binning" << fName << "' found... aborting." << endl;
705 return kFALSE;
706 }
707 }
708 if (binsx->HasTitle())
709 fHist->SetXTitle(binsx->GetTitle());
710 if (binsx->IsLogarithmic())
711 fHist->SetBit(kIsLogx);
712 }
713 }
714
715 // PreProcess existing fData members
716 for (int i=0; i<4; i++)
717 if (fData[i] && !fData[i]->PreProcess(plist))
718 return kFALSE;
719
720 if (fWeight && !fWeight->PreProcess(plist))
721 return kFALSE;
722
723 TString title(fDimension>0?"Histogram":"Profile");
724 title += " for ";
725 title += name;
726 title += Form(" (%dD)", TMath::Abs(fDimension));
727
728 fHist->SetName(name);
729 fHist->SetTitle(fTitle==gsDefTitle ? title : fTitle);
730 fHist->SetDirectory(0);
731
732 // This is for the case we have set lables
733 const MBinning def(1, 0, 1);
734 if (!binsx)
735 binsx = const_cast<MBinning*>(&def);
736 if (!binsy)
737 binsy = const_cast<MBinning*>(&def);
738 if (!binsz)
739 binsz = const_cast<MBinning*>(&def);
740
741 // set binning
742 switch (TMath::Abs(fDimension))
743 {
744 case 1:
745 SetBinning(*fHist, *binsx);
746 return kTRUE;
747 case 2:
748 SetBinning(static_cast<TH2&>(*fHist), *binsx, *binsy);
749 return kTRUE;
750 case 3:
751 SetBinning(static_cast<TH3&>(*fHist), *binsx, *binsy, *binsz);
752 return kTRUE;
753 }
754
755 *fLog << err << "ERROR - MH3 has " << TMath::Abs(fDimension) << " dimensions!" << endl;
756 return kFALSE;
757}
758
759// --------------------------------------------------------------------------
760//
761// Set the name of the histogram ant the MH3 container
762//
763void MH3::SetName(const char *name)
764{
765 if (fHist)
766 {
767 if (gPad)
768 {
769 const TString pfx(MString::Format("%sProfX", fHist->GetName()));
770 const TString pfy(MString::Format("%sProfY", fHist->GetName()));
771
772 TProfile *p = 0;
773 if ((p=dynamic_cast<TProfile*>(gPad->FindObject(pfx))))
774 p->SetName(MString::Format("%sProfX", name));
775 if ((p=dynamic_cast<TProfile*>(gPad->FindObject(pfy))))
776 p->SetName(MString::Format("%sProfY", name));
777 }
778
779 fHist->SetName(name);
780 fHist->SetDirectory(0);
781
782 }
783 MParContainer::SetName(name);
784}
785
786// --------------------------------------------------------------------------
787//
788// Set the title of the histogram ant the MH3 container
789//
790void MH3::SetTitle(const char *title)
791{
792 if (fHist)
793 fHist->SetTitle(title);
794 MParContainer::SetTitle(title);
795}
796
797// --------------------------------------------------------------------------
798//
799// Fills the one, two or three data members into our histogram
800//
801Int_t MH3::Fill(const MParContainer *par, const Stat_t ww)
802{
803 // Get Information about labels (UInt_t, to supress warning about
804 // unhandeled cases in switch)
805 const UInt_t type = GetLabels();
806
807 // Get values for axis
808 Double_t x=0;
809 Double_t y=0;
810 Double_t z=0;
811 Double_t t=0;
812 Double_t w=ww;
813
814 switch (fDimension)
815 {
816 case -3:
817 t = fData[3]->GetValue()*fScale[3];
818 case -2:
819 case 3:
820 z = fData[2]->GetValue()*fScale[2];
821 case -1:
822 case 2:
823 y = fData[1]->GetValue()*fScale[1];
824 case 1:
825 x = fData[0]->GetValue()*fScale[0];
826 }
827
828 if (fWeight)
829 w *= fWeight->GetValue();
830
831 // If label option is set, convert value to label
832 TString labelx, labely, labelz;
833 if (type&kLabelsX)
834 labelx = GetLabel(0, x);
835 if (type&kLabelsY)
836 labely = GetLabel(1, y);
837 if (type&kLabelsZ)
838 labelz = GetLabel(2, z);
839
840 // Fill histogram
841 switch (fDimension)
842 {
843 case 3:
844 switch (type)
845 {
846 case kNoLabels:
847 static_cast<TH3*>(fHist)->Fill(x, y, z, w);
848 return kTRUE;
849 case kLabelsX:
850#if ROOT_VERSION_CODE < ROOT_VERSION(5,32,01)
851 static_cast<TH3*>(fHist)->Fill(labelx, y, z);
852#endif
853 return kTRUE;
854 case kLabelsY:
855 static_cast<TH3*>(fHist)->Fill(x, labely, z, w);
856 return kTRUE;
857 case kLabelsZ:
858 static_cast<TH3*>(fHist)->Fill(x, y, labelz, w);
859 return kTRUE;
860 case kLabelsXY:
861 static_cast<TH3*>(fHist)->Fill(labelx, labely, z, w);
862 return kTRUE;
863 case kLabelsXZ:
864 static_cast<TH3*>(fHist)->Fill(labelx, y, labelz, w);
865 return kTRUE;
866 case kLabelsYZ:
867 static_cast<TH3*>(fHist)->Fill(x, labely, labelz, w);
868 return kTRUE;
869 case kLabelsXYZ:
870 static_cast<TH3*>(fHist)->Fill(labelx, labely, labelz, w);
871 return kTRUE;
872 }
873 break;
874 case 2:
875 switch (type)
876 {
877 case kNoLabels:
878 static_cast<TH2*>(fHist)->Fill(x, y, w);
879 return kTRUE;
880 case kLabelsX:
881 static_cast<TH2*>(fHist)->Fill(x, labely, w);
882 return kTRUE;
883 case kLabelsY:
884 static_cast<TH2*>(fHist)->Fill(labelx, y, w);
885 return kTRUE;
886 case kLabelsXY:
887 static_cast<TH2*>(fHist)->Fill(labelx, labely, w);
888 return kTRUE;
889 }
890 break;
891 case 1:
892 switch (type)
893 {
894 case kNoLabels:
895 fHist->Fill(x, w);
896 return kTRUE;
897 case kLabelsX:
898 fHist->Fill(labelx, w);
899 return kTRUE;
900 }
901 break;
902 case -1:
903 switch (type)
904 {
905 case kNoLabels:
906 static_cast<TProfile*>(fHist)->Fill(x, y, w);
907 return kTRUE;
908 case kLabelsX:
909 static_cast<TProfile*>(fHist)->Fill(labelx, y, w);
910 return kTRUE;
911 }
912 break;
913 case -2:
914 switch (type)
915 {
916 case kNoLabels:
917 static_cast<TProfile2D*>(fHist)->Fill(x, y, z, w);
918 return kTRUE;
919 case kLabelsX:
920 static_cast<TProfile2D*>(fHist)->Fill(labelx, y, z, 1.);
921 return kTRUE;
922 case kLabelsY:
923 static_cast<TProfile2D*>(fHist)->Fill(x, labely, z, 1.);
924 return kTRUE;
925 case kLabelsXY:
926 static_cast<TProfile2D*>(fHist)->Fill(labelx, labely, z, 1.);
927 return kTRUE;
928 }
929 break;
930 case -3:
931 switch (type)
932 {
933 case kNoLabels:
934 static_cast<TProfile3D*>(fHist)->Fill(x, y, z, t, w);
935 return kTRUE;
936 default:
937 *fLog << err << "ERROR - Labels not supported in TProfile3D." << endl;
938 }
939 break;
940 }
941
942 *fLog << err << "MH3::Fill: ERROR - A fatal error occured." << endl;
943 return kERROR;
944}
945
946// --------------------------------------------------------------------------
947//
948// If an auto range bit is set the histogram range of the corresponding
949// axis is set to show only the non-empty bins (with a single empty bin
950// on both sides)
951//
952Bool_t MH3::Finalize()
953{
954 DeflateLabels();
955
956 Bool_t autorangex=TESTBIT(fStyleBits, 0);
957 Bool_t autorangey=TESTBIT(fStyleBits, 1);
958 //Bool_t autorangez=TESTBIT(fStyleBits, 2);
959
960 Int_t lo, hi;
961 if (autorangex)
962 {
963 GetRangeX(*fHist, lo, hi);
964 fHist->GetXaxis()->SetRange(lo-2, hi+1);
965 }
966 if (autorangey)
967 {
968 GetRangeY(*fHist, lo, hi);
969 fHist->GetYaxis()->SetRange(lo-2, hi+1);
970 }
971 /*
972 if (autorangez)
973 {
974 GetRangeZ(*fHist, lo, hi);
975 fHist->GetZaxis()->SetRange(lo-2, hi+1);
976 }
977 */
978 return kTRUE;
979}
980
981// --------------------------------------------------------------------------
982//
983// Apply the conversion function to the contents stored in fHist and
984// store the result in h.
985//
986// In the case of a TProfile we keep it a TProfile to keep the mean and
987// rms in y displayed. To get the error correctly we have to reverse
988// the calculation done in TProfile::GetBinError of course.
989//
990void MH3::Convert(TH1 &h) const
991{
992 const Bool_t prof = h.InheritsFrom(TProfile::Class()) || h.InheritsFrom(TProfile2D::Class()) || h.InheritsFrom(TProfile3D::Class());
993
994 for (Int_t z=0; z<=h.GetNbinsZ()+1; z++)
995 for (Int_t y=0; y<=h.GetNbinsY()+1; y++)
996 for (Int_t x=0; x<=h.GetNbinsX()+1; x++)
997 {
998 h.SetBinContent(x, y, z, fConversion->Eval(fHist->GetBinContent(x, y, z)));
999
1000 if (prof)
1001 h.SetBinError(x, y, z, TMath::Hypot(fConversion->Eval(fHist->GetBinContent(x, y, z)),
1002 fConversion->Eval(fHist->GetBinError( x, y, z))));
1003 else
1004 h.SetBinError(x, y, z, fConversion->Eval(fHist->GetBinError(x, y, z)));
1005 }
1006
1007 TProfile *p1 = dynamic_cast<TProfile*>(fHist);
1008 if (p1)
1009 for (Int_t i=0; i<p1->GetSize(); i++)
1010 static_cast<TProfile&>(h).SetBinEntries(i, p1->GetBinEntries(i)>0 ? 1 : 0);
1011
1012 TProfile2D *p2 = dynamic_cast<TProfile2D*>(fHist);
1013 if (p2)
1014 for (Int_t i=0; i<p2->GetSize(); i++)
1015 static_cast<TProfile2D&>(h).SetBinEntries(i, p2->GetBinEntries(i)>0 ? 1 : 0);
1016
1017 TProfile3D *p3 = dynamic_cast<TProfile3D*>(fHist);
1018 if (p3)
1019 for (Int_t i=0; i<p3->GetSize(); i++)
1020 static_cast<TProfile3D&>(h).SetBinEntries(i, p3->GetBinEntries(i)>0 ? 1 : 0);
1021}
1022
1023// --------------------------------------------------------------------------
1024//
1025// FIXME
1026//
1027void MH3::Paint(Option_t *o)
1028{
1029 TProfile *p=0;
1030
1031 if (TMath::Abs(fDimension)==2)
1032 MH::SetPalette("pretty");
1033
1034 if (fConversion)
1035 {
1036 TH1 *h = 0;
1037 if ((h=dynamic_cast<TH1*>(gPad->FindObject(fHist->GetName()))))
1038 Convert(*h);
1039 }
1040
1041 const TString pfx(MString::Format("%sProfX", fHist->GetName()));
1042 if ((p=dynamic_cast<TProfile*>(gPad->FindObject(pfx))))
1043 {
1044 Int_t col = p->GetLineColor();
1045 p = ((TH2*)fHist)->ProfileX(pfx, -1, -1, "s");
1046 p->SetLineColor(col);
1047 }
1048
1049 const TString pfy(MString::Format("%sProfY", fHist->GetName()));
1050 if ((p=dynamic_cast<TProfile*>(gPad->FindObject(pfy))))
1051 {
1052 Int_t col = p->GetLineColor();
1053 p = ((TH2*)fHist)->ProfileY(pfy, -1, -1, "s");
1054 p->SetLineColor(col);
1055 }
1056/*
1057 if (fHist->TestBit(kIsLogx) && fHist->GetEntries()>0) gPad->SetLogx();
1058 if (fHist->TestBit(kIsLogy) && fHist->GetEntries()>0) gPad->SetLogy();
1059 if (fHist->TestBit(kIsLogz) && fHist->GetEntries()>0) gPad->SetLogz();
1060 */
1061}
1062
1063// --------------------------------------------------------------------------
1064//
1065// If Xmax is < 3000*Xmin SetMoreLogLabels is called. If Xmax<5000
1066// the exponent is switched off (SetNoExponent)
1067//
1068void MH3::HandleLogAxis(TAxis &axe) const
1069{
1070 if (axe.GetXmax()>3000*axe.GetXmin())
1071 return;
1072
1073 axe.SetMoreLogLabels();
1074 if (axe.GetXmax()<5000)
1075 axe.SetNoExponent();
1076}
1077
1078// --------------------------------------------------------------------------
1079//
1080// Creates a new canvas and draws the histogram into it.
1081//
1082// Possible options are:
1083// PROFX: Draw a x-profile into the histogram (for 2D histograms only)
1084// PROFY: Draw a y-profile into the histogram (for 2D histograms only)
1085// ONLY: Draw the profile histogram only (for 2D histograms only)
1086// BLUE: Draw the profile in blue color instead of the histograms
1087// line color
1088//
1089// If the kIsLog?-Bit is set the axis is displayed lkogarithmically.
1090// eg this is set when applying a logarithmic MBinning
1091//
1092// Be careful: The histogram belongs to this object and won't get deleted
1093// together with the canvas.
1094//
1095void MH3::Draw(Option_t *opt)
1096{
1097 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(fHist);
1098 pad->SetBorderMode(0);
1099 pad->SetGridx();
1100 pad->SetGridy();
1101
1102 AppendPad();
1103
1104 if (fHist->TestBit(kIsLogx))
1105 {
1106 pad->SetLogx();
1107 HandleLogAxis(*fHist->GetXaxis());
1108 }
1109 if (fHist->TestBit(kIsLogy))
1110 {
1111 pad->SetLogy();
1112 HandleLogAxis(*fHist->GetYaxis());
1113 }
1114 if (fHist->TestBit(kIsLogz))
1115 {
1116 pad->SetLogz();
1117 HandleLogAxis(*fHist->GetZaxis());
1118 }
1119
1120 fHist->SetFillStyle(4000);
1121
1122 TString str(opt);
1123 str.ToLower();
1124 str.ReplaceAll(" ", "");
1125
1126 if (GetLabels())
1127 {
1128 if (str.IsNull() && fDimension==2)
1129 str = "colz";
1130
1131 if (str.Contains("box", TString::kIgnoreCase) && fDimension==2)
1132 fHist->SetLineColor(kBlue);
1133 }
1134
1135 const Bool_t only = str.Contains("only") && TMath::Abs(fDimension)==2;
1136 const Bool_t same = str.Contains("same") && TMath::Abs(fDimension)<3;
1137 const Bool_t blue = str.Contains("blue") && TMath::Abs(fDimension)==2;
1138 const Bool_t profx = str.Contains("profx") && TMath::Abs(fDimension)==2;
1139 const Bool_t profy = str.Contains("profy") && TMath::Abs(fDimension)==2;
1140
1141 str.ReplaceAll("only", "");
1142 str.ReplaceAll("blue", "");
1143 str.ReplaceAll("profx", "");
1144 str.ReplaceAll("profy", "");
1145
1146 if (same && TMath::Abs(fDimension)==1)
1147 {
1148 fHist->SetLineColor(kBlue);
1149 fHist->SetMarkerColor(kBlue);
1150 }
1151
1152
1153 TH1 *h = fHist;
1154
1155 if (fConversion)
1156 {
1157 h = static_cast<TH1*>(fHist->Clone());
1158
1159 h->SetDirectory(0);
1160 h->SetBit(kCanDelete);
1161
1162 Convert(*h);
1163 }
1164
1165 // FIXME: We may have to remove all our own options from str!
1166 if (!only)
1167 h->Draw(str);
1168
1169 TProfile *p=0;
1170 if (profx)
1171 {
1172 const TString pfx(MString::Format("%sProfX", h->GetName()));
1173
1174 if (same && (p=dynamic_cast<TProfile*>(gPad->FindObject(pfx))))
1175 *fLog << warn << "TProfile " << pfx << " already in pad." << endl;
1176
1177 p = ((TH2*)h)->ProfileX(pfx, -1, -1, "s");
1178 p->UseCurrentStyle();
1179 p->SetLineColor(blue ? kBlue : h->GetLineColor());
1180 p->SetBit(kCanDelete);
1181 p->SetDirectory(NULL);
1182 p->SetXTitle(h->GetXaxis()->GetTitle());
1183 p->SetYTitle(h->GetYaxis()->GetTitle());
1184 p->Draw(only&&!same?"":"same");
1185 }
1186 if (profy)
1187 {
1188 const TString pfy(MString::Format("%sProfY", h->GetName()));
1189
1190 if (same && (p=dynamic_cast<TProfile*>(gPad->FindObject(pfy))))
1191 *fLog << warn << "TProfile " << pfy << " already in pad." << endl;
1192
1193 p = ((TH2*)h)->ProfileY(pfy, -1, -1, "s");
1194 p->UseCurrentStyle();
1195 p->SetLineColor(blue ? kBlue : h->GetLineColor());
1196 p->SetBit(kCanDelete);
1197 p->SetDirectory(NULL);
1198 p->SetYTitle(h->GetXaxis()->GetTitle());
1199 p->SetXTitle(h->GetYaxis()->GetTitle());
1200 p->Draw(only&&!same?"":"same");
1201 }
1202
1203 //AppendPad("log");
1204}
1205
1206// --------------------------------------------------------------------------
1207//
1208// Implementation of SavePrimitive. Used to write the call to a constructor
1209// to a macro. In the original root implementation it is used to write
1210// gui elements to a macro-file.
1211//
1212void MH3::StreamPrimitive(ostream &out) const
1213{
1214 TString name = GetUniqueName();
1215
1216 out << " MH3 " << name << "(\"";
1217 out << fData[0]->GetRule() << "\"";
1218 if (fDimension>1 || fDimension<0)
1219 out << ", \"" << fData[1]->GetRule() << "\"";
1220 if (fDimension>2 || fDimension<-1)
1221 out << ", \"" << fData[2]->GetRule() << "\"";
1222
1223 out << ");" << endl;
1224
1225 if (fName!=gsDefName)
1226 out << " " << name << ".SetName(\"" << fName << "\");" << endl;
1227
1228 if (fTitle!=gsDefTitle)
1229 out << " " << name << ".SetTitle(\"" << fTitle << "\");" << endl;
1230
1231 if (fWeight)
1232 out << " " << name << ".SetWeight(\"" << fWeight->GetRule() << "\");" << endl;
1233
1234 switch (fDimension)
1235 {
1236 case -3:
1237 if (fScale[3]!=1)
1238 out << " " << name << ".SetScaleT(" << fScale[3] << ");" << endl;
1239 case -2:
1240 case 3:
1241 if (fScale[2]!=1)
1242 out << " " << name << ".SetScaleZ(" << fScale[2] << ");" << endl;
1243 case -1:
1244 case 2:
1245 if (fScale[1]!=1)
1246 out << " " << name << ".SetScaleY(" << fScale[1] << ");" << endl;
1247 case 1:
1248 if (fScale[0]!=1)
1249 out << " " << name << ".SetScaleX(" << fScale[0] << ");" << endl;
1250 }
1251}
1252
1253MH3::Type_t MH3::GetType() const
1254{
1255 switch (fDimension)
1256 {
1257 case -1:
1258 return TString(static_cast<TProfile*>(fHist)->GetErrorOption())=="s" ? kProfileSpread : kProfile;
1259 case -2:
1260 return TString(static_cast<TProfile2D*>(fHist)->GetErrorOption())=="s" ? kProfileSpread : kProfile;
1261 case -3:
1262 return TString(static_cast<TProfile3D*>(fHist)->GetErrorOption())=="s" ? kProfileSpread : kProfile;
1263 }
1264 return kHistogram;
1265}
1266
1267// --------------------------------------------------------------------------
1268//
1269// Used to rebuild a MH3 object of the same type (data members,
1270// dimension, ...)
1271//
1272MParContainer *MH3::New() const
1273{
1274 // FIXME: TREAT THE NEW OPTIONS CORRECTLY (PROFILE, LABELS)
1275
1276 MH3 *h = NULL;
1277
1278 if (fData[0] == NULL)
1279 h=new MH3(fDimension);
1280 else
1281 switch (fDimension)
1282 {
1283 case 1:
1284 h=new MH3(fData[0]->GetRule());
1285 break;
1286 case 2:
1287 case -1:
1288 h=new MH3(fData[0]->GetRule(), fData[1]->GetRule(), GetType());
1289 break;
1290 case 3:
1291 case -2:
1292 h=new MH3(fData[0]->GetRule(), fData[1]->GetRule(), fData[2]->GetRule(), GetType());
1293 break;
1294 case -3:
1295 h=new MH3(fData[0]->GetRule(), fData[1]->GetRule(), fData[2]->GetRule(), fData[3]->GetRule());
1296 break;
1297 }
1298
1299 switch (fDimension)
1300 {
1301 case -3:
1302 h->SetScaleZ(fScale[3]);
1303 case -2:
1304 case 3:
1305 h->SetScaleZ(fScale[2]);
1306 case 2:
1307 case -1:
1308 h->SetScaleY(fScale[1]);
1309 case 1:
1310 h->SetScaleX(fScale[0]);
1311 }
1312
1313 if (fWeight)
1314 h->SetWeight(fWeight->GetRule());
1315
1316 return h;
1317}
1318
1319// --------------------------------------------------------------------------
1320//
1321// FIXME
1322//
1323TString MH3::GetRule(const Char_t axis) const
1324{
1325 switch (tolower(axis))
1326 {
1327 case 'x':
1328 case 'X':
1329 return fData[0] ? fData[0]->GetRule() : TString("");
1330 case 'y':
1331 case 'Y':
1332 return fData[1] ? fData[1]->GetRule() : TString("");
1333 case 'z':
1334 case 'Z':
1335 return fData[2] ? fData[2]->GetRule() : TString("");
1336 case 't':
1337 case 'T':
1338 return fData[3] ? fData[3]->GetRule() : TString("");
1339 case 'w':
1340 case 'W':
1341 return fWeight ? fWeight->GetRule() : TString("");
1342 default:
1343 return "<n/a>";
1344 }
1345}
1346
1347// --------------------------------------------------------------------------
1348//
1349// Returns the total number of bins in a histogram (excluding under- and
1350// overflow bins)
1351//
1352Int_t MH3::GetNbins() const
1353{
1354 Int_t num = 1;
1355
1356 switch (TMath::Abs(fDimension))
1357 {
1358 case 3:
1359 num *= fHist->GetNbinsZ()+2;
1360 case 2:
1361 num *= fHist->GetNbinsY()+2;
1362 case 1:
1363 num *= fHist->GetNbinsX()+2;
1364 }
1365
1366 return num;
1367}
1368
1369// --------------------------------------------------------------------------
1370//
1371// Returns the total number of bins in a histogram (excluding under- and
1372// overflow bins) Return -1 if bin is underflow or overflow
1373//
1374Int_t MH3::FindFixBin(Double_t x, Double_t y, Double_t z) const
1375{
1376 const TAxis &axex = *fHist->GetXaxis();
1377 const TAxis &axey = *fHist->GetYaxis();
1378 const TAxis &axez = *fHist->GetZaxis();
1379
1380 Int_t binz = 0;
1381 Int_t biny = 0;
1382 Int_t binx = 0;
1383
1384 switch (fDimension)
1385 {
1386 case 3:
1387 case -2:
1388 binz = axez.FindFixBin(z);
1389 if (binz>axez.GetLast() || binz<axez.GetFirst())
1390 return -1;
1391 case 2:
1392 case -1:
1393 biny = axey.FindFixBin(y);
1394 if (biny>axey.GetLast() || biny<axey.GetFirst())
1395 return -1;
1396 case 1:
1397 binx = axex.FindFixBin(x);
1398 if (binx<axex.GetFirst() || binx>axex.GetLast())
1399 return -1;
1400 }
1401
1402 const Int_t nx = fHist->GetNbinsX()+2;
1403 const Int_t ny = fHist->GetNbinsY()+2;
1404
1405 return binx + nx*(biny +ny*binz);
1406}
1407
1408// --------------------------------------------------------------------------
1409//
1410// Return the MObjLookup corresponding to the axis/character.
1411// Note that only lower-case charecters (x, y, z) are supported.
1412// If for the axis no labels were set, the corresponding
1413// InitLabels is called.
1414//
1415MObjLookup *MH3::GetLabels(char axe)
1416{
1417 if (!fHist)
1418 return 0;
1419
1420 TAxis *x = 0;
1421
1422 switch (axe)
1423 {
1424 case 'x':
1425 x = fHist->GetXaxis();
1426 break;
1427 case 'y':
1428 x = fHist->GetYaxis();
1429 break;
1430 case 'z':
1431 x = fHist->GetZaxis();
1432 break;
1433 }
1434
1435 if (!x)
1436 return 0;
1437
1438 const Int_t idx = axe-'x';
1439
1440 if (!x->GetLabels())
1441 switch (idx)
1442 {
1443 case 0:
1444 InitLabels(kLabelsX);
1445 break;
1446 case 1:
1447 InitLabels(kLabelsY);
1448 break;
1449 case 2:
1450 InitLabels(kLabelsZ);
1451 break;
1452 }
1453
1454 return &fLabels[idx];
1455}
1456
1457// --------------------------------------------------------------------------
1458//
1459// Set a default label which is used if no other is found in the list
1460// of labels. if a default was set already it is overwritten. If the
1461// axis has not yet been initialized to use labels it it now.
1462//
1463void MH3::DefaultLabel(char axe, const char *name)
1464{
1465 MObjLookup *arr = GetLabels(axe);
1466 if (!arr)
1467 return;
1468
1469 if (arr->GetDefault())
1470 {
1471 delete arr->GetDefault();
1472 arr->SetDefault(0);
1473 }
1474
1475 if (name)
1476 arr->SetDefault(new TObjString(name));
1477}
1478
1479// --------------------------------------------------------------------------
1480//
1481// Define a name for a label. More than one label can have the same
1482// name. If the axis has not yet been initialized to use labels
1483// it it now.
1484//
1485void MH3::DefineLabel(char axe, Int_t label, const char *name)
1486{
1487 MObjLookup *arr = GetLabels(axe);
1488
1489 if (!arr || !name)
1490 return;
1491
1492 if (arr->GetObj(label)!=arr->GetDefault())
1493 return;
1494
1495 arr->Add(label, name);
1496}
1497
1498// --------------------------------------------------------------------------
1499//
1500// Define names for labels, like
1501// 1=Trig;2=Cal;4=Ped;8=Lvl2
1502// More than one label can have the same name. If the axis has not
1503// yet been initialized to use labels it it now.
1504//
1505// A default cannot be set here. Use DefaultLabel instead.
1506//
1507void MH3::DefineLabels(char axe, const TString &labels)
1508{
1509 TObjArray *arr = labels.Tokenize(';');
1510
1511 for (int i=0; i<arr->GetEntries(); i++)
1512 {
1513 const char *s = (*arr)[0]->GetName();
1514 const char *v = strchr(s, '=');
1515
1516 if (v)
1517 DefineLabel(axe, atoi(s), v+1);
1518 }
1519
1520 delete arr;
1521}
1522
1523void MH3::RecursiveRemove(TObject *obj)
1524{
1525 if (obj==fHist)
1526 fHist = 0;
1527
1528 MH::RecursiveRemove(obj);
1529}
Note: See TracBrowser for help on using the repository browser.