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

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