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

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