source: trunk/MagicSoft/Mars/mhbase/MH3.cc@ 8813

Last change on this file since 8813 was 8709, checked in by tbretz, 17 years ago
*** empty log message ***
File size: 22.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-2007
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 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//
61// Axis titles
62// ===========
63//
64// 1) If no other title is given the rule for this axis is used.
65// 2) If the MBinning used for this axis has a non-default Title
66// (MBinning::HasTitle) this title is used for the corresponding axis
67// 3) If the MH3 has a non-default title (MH3::SetTitle called)
68// this title is set as the histogram title. It can be used to overwrite
69// the axis titles. For more information see TH1::SetTitle, eg.
70// SetTitle("MyHist;x[mm];y[cm];Counts");
71//
72//
73// For example:
74// MH3 myhist("MHillas.fLength");
75// myhist.SetName("MyHist");
76// myhist.SetScaleX(geomcam.GetConvMm2Deg()); //convert length to degree
77// MBinning bins("BinningMyHistX");
78// bins.SetEdges(10, 0, 150);
79// plist.AddToList(&myhist);
80// plist.AddToList(&bins);
81//
82//
83// Class Version 2:
84// ----------------
85// - MDataChain *fData[3]; // Object from which the data is filled
86// + MData *fData[3]; // Object from which the data is filled
87//
88// Class Version 3:
89// ----------------
90// - Byte_t fStyleBits
91// + MBinning fBins[3]
92//
93/////////////////////////////////////////////////////////////////////////////
94#include "MH3.h"
95
96#include <ctype.h> // tolower
97#include <fstream>
98
99//#include <TPad.h>
100#include <TStyle.h>
101#include <TCanvas.h>
102
103#include <TH2.h>
104#include <TH3.h>
105#include <TProfile.h>
106#include <TProfile2D.h>
107
108#include "MLog.h"
109#include "MLogManip.h"
110
111#include "MParList.h"
112#include "MBinning.h"
113#include "MDataPhrase.h"
114
115ClassImp(MH3);
116
117using namespace std;
118
119const TString MH3::gsDefName = "MH3";
120const TString MH3::gsDefTitle = "Container for a n-D Mars Histogram";
121
122// --------------------------------------------------------------------------
123//
124// Default constructor.
125//
126MH3::MH3(const unsigned int dim)
127 : fDimension(dim>3?3:dim), fHist(NULL), fStyleBits(0)
128{
129 switch (fDimension)
130 {
131 case 1:
132 fHist = new TH1D;
133 fHist->SetYTitle("Counts");
134 break;
135 case 2:
136 fHist = new TH2D;
137 fHist->SetZTitle("Counts");
138 break;
139 case 3:
140 fHist = new TH3D;
141 break;
142 }
143
144 fData[0] = NULL;
145 fData[1] = NULL;
146 fData[2] = NULL;
147
148 fBins[0] = NULL;
149 fBins[1] = NULL;
150 fBins[2] = NULL;
151
152 fName = gsDefName;
153 fTitle = gsDefTitle;
154
155 if (fHist)
156 {
157 fHist->SetDirectory(NULL);
158 fHist->UseCurrentStyle();
159 }
160
161 fScale[0] = 1;
162 fScale[1] = 1;
163 fScale[2] = 1;
164}
165
166// --------------------------------------------------------------------------
167//
168// Creates an TH1D. memberx is filled into the X-bins. For a more detailed
169// description see the class description above.
170//
171MH3::MH3(const char *memberx)
172 : fDimension(1), fStyleBits(0)
173{
174 fHist = new TH1D;
175
176 fData[0] = new MDataPhrase(memberx);
177 fData[1] = NULL;
178 fData[2] = NULL;
179
180 fBins[0] = NULL;
181 fBins[1] = NULL;
182 fBins[2] = NULL;
183
184 fName = gsDefName;
185 fTitle = gsDefTitle;
186
187 fHist->UseCurrentStyle();
188 fHist->SetDirectory(NULL);
189 fHist->SetYTitle("Counts");
190
191 fScale[0] = 1;
192 fScale[1] = 1;
193 fScale[2] = 1;
194}
195
196MH3::MH3(const TH1 &h1)
197 : fDimension(1), fStyleBits(0)
198{
199 if (h1.InheritsFrom(TH3::Class()))
200 fDimension = 3;
201 if (h1.InheritsFrom(TH2::Class()))
202 fDimension = 2;
203
204 fData[0] = NULL;
205 fData[1] = NULL;
206 fData[2] = NULL;
207
208 fBins[0] = NULL;
209 fBins[1] = NULL;
210 fBins[2] = NULL;
211
212 switch (fDimension)
213 {
214 case 3:
215 fData[2] = new MDataPhrase(h1.GetZaxis()->GetTitle());
216 case 2:
217 fData[1] = new MDataPhrase(h1.GetYaxis()->GetTitle());
218 case 1:
219 fData[0] = new MDataPhrase(h1.GetXaxis()->GetTitle());
220 }
221
222 fName = gsDefName;
223 fTitle = gsDefTitle;
224
225 fHist = (TH1*)h1.Clone();
226 fHist->SetDirectory(NULL);
227
228 fScale[0] = 1;
229 fScale[1] = 1;
230 fScale[2] = 1;
231}
232
233// --------------------------------------------------------------------------
234//
235// Creates an TH2D. memberx is filled into the X-bins. membery is filled
236// into the Y-bins. For a more detailed description see the class
237// description above.
238//
239MH3::MH3(const char *memberx, const char *membery)
240 : fDimension(2), fStyleBits(0)
241{
242 fHist = new TH2D;
243
244 fData[0] = new MDataPhrase(memberx);
245 fData[1] = new MDataPhrase(membery);
246 fData[2] = NULL;
247
248 fBins[0] = NULL;
249 fBins[1] = NULL;
250 fBins[2] = NULL;
251
252 fName = gsDefName;
253 fTitle = gsDefTitle;
254
255 fHist->UseCurrentStyle();
256 fHist->SetDirectory(NULL);
257 fHist->SetZTitle("Counts");
258
259 fScale[0] = 1;
260 fScale[1] = 1;
261 fScale[2] = 1;
262}
263
264// --------------------------------------------------------------------------
265//
266// Creates an TH3D. memberx is filled into the X-bins. membery is filled
267// into the Y-bins. membery is filled into the Z-bins. For a more detailed
268// description see the class description above.
269//
270MH3::MH3(const char *memberx, const char *membery, const char *memberz)
271 : fDimension(3), fStyleBits(0)
272{
273 fHist = new TH3D;
274
275 fData[0] = new MDataPhrase(memberx);
276 fData[1] = new MDataPhrase(membery);
277 fData[2] = new MDataPhrase(memberz);
278
279 fBins[0] = NULL;
280 fBins[1] = NULL;
281 fBins[2] = NULL;
282
283 fName = gsDefName;
284 fTitle = gsDefTitle;
285
286 fHist->UseCurrentStyle();
287 fHist->SetDirectory(NULL);
288
289 fScale[0] = 1;
290 fScale[1] = 1;
291 fScale[2] = 1;
292}
293
294// --------------------------------------------------------------------------
295//
296// Deletes the histogram
297//
298MH3::~MH3()
299{
300 delete fHist;
301
302 for (int i=0; i<3; i++)
303 if (fData[i])
304 delete fData[i];
305}
306
307// --------------------------------------------------------------------------
308//
309// Return the data members used by the data chain to be used in
310// MTask::AddBranchToList
311//
312TString MH3::GetDataMember() const
313{
314 TString str=fData[0]->GetDataMember();
315 if (fData[1])
316 {
317 str += ";";
318 str += fData[1]->GetDataMember();
319 }
320 if (fData[2])
321 {
322 str += ";";
323 str += fData[2]->GetDataMember();
324 }
325 return str;
326}
327
328// --------------------------------------------------------------------------
329//
330// Setup the Binning for the histograms automatically if the correct
331// instances of MBinning are found in the parameter list
332// For a more detailed description see class description above.
333//
334Bool_t MH3::SetupFill(const MParList *plist)
335{
336 // reset histogram (necessary if the same eventloop is run more than once)
337 fHist->Reset();
338
339 // Tokenize name into name and binnings names
340 TObjArray *tok = fName.Tokenize(";");
341
342 const TString name = (*tok)[0] ? (*tok)[0]->GetName() : fName.Data();
343
344 TString bx = (*tok)[1] ? (*tok)[1]->GetName() : Form("%sX", name.Data());
345 TString by = (*tok)[2] ? (*tok)[2]->GetName() : Form("%sY", name.Data());
346 TString bz = (*tok)[3] ? (*tok)[3]->GetName() : Form("%sZ", name.Data());
347
348 bx.Prepend("Binning");
349 by.Prepend("Binning");
350 bz.Prepend("Binning");
351
352 delete tok;
353
354
355 MBinning *binsx = NULL;
356 MBinning *binsy = NULL;
357 MBinning *binsz = NULL;
358
359 switch (fDimension)
360 {
361 case 3:
362 binsz = fBins[2] ? fBins[2] : (MBinning*)plist->FindObject(bz, "MBinning");
363 if (!binsz)
364 {
365 *fLog << err << dbginf << "MBinning '" << bz << "' not found... aborting." << endl;
366 return kFALSE;
367 }
368 if (fData[2] && !fData[2]->PreProcess(plist))
369 return kFALSE;
370 if (fData[2])
371 fHist->SetZTitle(fData[2]->GetTitle());
372 if (binsz->HasTitle())
373 fHist->SetZTitle(binsz->GetTitle());
374 if (binsz->IsLogarithmic())
375 fHist->SetBit(kIsLogz);
376 case 2:
377 binsy = fBins[1] ? fBins[1] : (MBinning*)plist->FindObject(by, "MBinning");
378 if (!binsy)
379 {
380 *fLog << err << dbginf << "MBinning '" << by << "' not found... aborting." << endl;
381 return kFALSE;
382 }
383 if (fData[1] && !fData[1]->PreProcess(plist))
384 return kFALSE;
385 if (fData[1])
386 fHist->SetYTitle(fData[1]->GetTitle());
387 if (binsy->HasTitle())
388 fHist->SetYTitle(binsy->GetTitle());
389 if (binsy->IsLogarithmic())
390 fHist->SetBit(kIsLogy);
391 case 1:
392 binsx = fBins[0] ? fBins[0] : (MBinning*)plist->FindObject(bx, "MBinning");
393 if (!binsx)
394 {
395 if (fDimension==1)
396 binsx = (MBinning*)plist->FindObject("Binning"+fName, "MBinning");
397
398 if (!binsx)
399 {
400 *fLog << err << dbginf << "Neither '" << bx << "' nor '" << binsx << fName << "' found... aborting." << endl;
401 return kFALSE;
402 }
403 }
404 if (fData[0] && !fData[0]->PreProcess(plist))
405 return kFALSE;
406 if (fData[0]!=NULL)
407 fHist->SetXTitle(fData[0]->GetTitle());
408 if (binsx->HasTitle())
409 fHist->SetXTitle(binsx->GetTitle());
410 if (binsx->IsLogarithmic())
411 fHist->SetBit(kIsLogx);
412 }
413
414 TString title("Histogram for ");
415 title += name;
416 title += Form(" (%dD)", fDimension);
417
418 fHist->SetName(name);
419 fHist->SetTitle(fTitle==gsDefTitle ? title : fTitle);
420 fHist->SetDirectory(0);
421
422 switch (fDimension)
423 {
424 case 1:
425 SetBinning(fHist, binsx);
426 return kTRUE;
427 case 2:
428 SetBinning((TH2*)fHist, binsx, binsy);
429 return kTRUE;
430 case 3:
431 SetBinning((TH3*)fHist, binsx, binsy, binsz);
432 return kTRUE;
433 }
434
435 *fLog << err << "ERROR - MH3 has " << fDimension << " dimensions!" << endl;
436 return kFALSE;
437}
438
439// --------------------------------------------------------------------------
440//
441// Set the name of the histogram ant the MH3 container
442//
443void MH3::SetName(const char *name)
444{
445 if (fHist)
446 {
447 if (gPad)
448 {
449 const TString pfx(Form("%sProfX", fHist->GetName()));
450 const TString pfy(Form("%sProfY", fHist->GetName()));
451
452 TProfile *p = 0;
453 if ((p=dynamic_cast<TProfile*>(gPad->FindObject(pfx))))
454 p->SetName(Form("%sProfX", name));
455 if ((p=dynamic_cast<TProfile*>(gPad->FindObject(pfy))))
456 p->SetName(Form("%sProfY", name));
457 }
458
459 fHist->SetName(name);
460 fHist->SetDirectory(0);
461
462 }
463 MParContainer::SetName(name);
464}
465
466// --------------------------------------------------------------------------
467//
468// Set the title of the histogram ant the MH3 container
469//
470void MH3::SetTitle(const char *title)
471{
472 if (fHist)
473 fHist->SetTitle(title);
474 MParContainer::SetTitle(title);
475}
476
477// --------------------------------------------------------------------------
478//
479// Fills the one, two or three data members into our histogram
480//
481Bool_t MH3::Fill(const MParContainer *par, const Stat_t w)
482{
483 Double_t x=0;
484 Double_t y=0;
485 Double_t z=0;
486
487 switch (fDimension)
488 {
489 case 3:
490 z = fData[2]->GetValue()*fScale[2];
491 case 2:
492 y = fData[1]->GetValue()*fScale[1];
493 case 1:
494 x = fData[0]->GetValue()*fScale[0];
495 }
496
497 switch (fDimension)
498 {
499 case 3:
500 ((TH3*)fHist)->Fill(x, y, z, w);
501 return kTRUE;
502 case 2:
503 ((TH2*)fHist)->Fill(x, y, w);
504 return kTRUE;
505 case 1:
506 fHist->Fill(x, w);
507 return kTRUE;
508 }
509
510 return kFALSE;
511}
512
513// --------------------------------------------------------------------------
514//
515// If an auto range bit is set the histogram range of the corresponding
516// axis is set to show only the non-empty bins (with a single empty bin
517// on both sides)
518//
519Bool_t MH3::Finalize()
520{
521 Bool_t autorangex=TESTBIT(fStyleBits, 0);
522 Bool_t autorangey=TESTBIT(fStyleBits, 1);
523 //Bool_t autorangez=TESTBIT(fStyleBits, 2);
524
525 Int_t lo, hi;
526
527 if (autorangex)
528 {
529 GetRangeX(*fHist, lo, hi);
530 fHist->GetXaxis()->SetRange(lo-2, hi+1);
531 }
532 if (autorangey)
533 {
534 GetRangeY(*fHist, lo, hi);
535 fHist->GetYaxis()->SetRange(lo-2, hi+1);
536 }
537 /*
538 if (autorangez)
539 {
540 GetRangeZ(*fHist, lo, hi);
541 fHist->GetZaxis()->SetRange(lo-2, hi+1);
542 }
543 */
544 return kTRUE;
545}
546
547void MH3::Paint(Option_t *o)
548{
549 TProfile *p=0;
550
551 if (fDimension==2)
552 MH::SetPalette("pretty");
553
554 const TString pfx(Form("%sProfX", fHist->GetName()));
555 if ((p=dynamic_cast<TProfile*>(gPad->FindObject(pfx))))
556 {
557 Int_t col = p->GetLineColor();
558 p = ((TH2*)fHist)->ProfileX(pfx, -1, -1, "s");
559 p->SetLineColor(col);
560 }
561
562 const TString pfy(Form("%sProfY", fHist->GetName()));
563 if ((p=dynamic_cast<TProfile*>(gPad->FindObject(pfy))))
564 {
565 Int_t col = p->GetLineColor();
566 p = ((TH2*)fHist)->ProfileY(pfy, -1, -1, "s");
567 p->SetLineColor(col);
568 }
569/*
570 if (fHist->TestBit(kIsLogx) && fHist->GetEntries()>0) gPad->SetLogx();
571 if (fHist->TestBit(kIsLogy) && fHist->GetEntries()>0) gPad->SetLogy();
572 if (fHist->TestBit(kIsLogz) && fHist->GetEntries()>0) gPad->SetLogz();
573 */
574}
575
576void MH3::HandleLogAxis(TAxis &axe) const
577{
578 if (axe.GetXmax()>3000*axe.GetXmin())
579 return;
580
581 axe.SetMoreLogLabels();
582 if (axe.GetXmax()<5000)
583 axe.SetNoExponent();
584}
585
586// --------------------------------------------------------------------------
587//
588// Creates a new canvas and draws the histogram into it.
589//
590// Possible options are:
591// PROFX: Draw a x-profile into the histogram (for 2D histograms only)
592// PROFY: Draw a y-profile into the histogram (for 2D histograms only)
593// ONLY: Draw the profile histogram only (for 2D histograms only)
594// BLUE: Draw the profile in blue color instead of the histograms
595// line color
596//
597// If the kIsLog?-Bit is set the axis is displayed lkogarithmically.
598// eg this is set when applying a logarithmic MBinning
599//
600// Be careful: The histogram belongs to this object and won't get deleted
601// together with the canvas.
602//
603void MH3::Draw(Option_t *opt)
604{
605 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(fHist);
606 pad->SetBorderMode(0);
607 pad->SetGridx();
608 pad->SetGridy();
609
610 if (fHist->TestBit(kIsLogx))
611 {
612 pad->SetLogx();
613 HandleLogAxis(*fHist->GetXaxis());
614 }
615 if (fHist->TestBit(kIsLogy))
616 {
617 pad->SetLogy();
618 HandleLogAxis(*fHist->GetYaxis());
619 }
620 if (fHist->TestBit(kIsLogz))
621 {
622 pad->SetLogz();
623 HandleLogAxis(*fHist->GetZaxis());
624 }
625
626 fHist->SetFillStyle(4000);
627
628 TString str(opt);
629 str.ToLower();
630
631 const Bool_t only = str.Contains("only") && fDimension==2;
632 const Bool_t same = str.Contains("same") && fDimension<3;
633 const Bool_t blue = str.Contains("blue") && fDimension==2;
634 const Bool_t profx = str.Contains("profx") && fDimension==2;
635 const Bool_t profy = str.Contains("profy") && fDimension==2;
636
637 str.ReplaceAll("only", "");
638 str.ReplaceAll("blue", "");
639 str.ReplaceAll("profx", "");
640 str.ReplaceAll("profy", "");
641 str.ReplaceAll(" ", "");
642
643 if (same && fDimension==1)
644 {
645 fHist->SetLineColor(kBlue);
646 fHist->SetMarkerColor(kBlue);
647 }
648
649 // FIXME: We may have to remove all our own options from str!
650 if (!only)
651 fHist->Draw(str);
652
653 AppendPad();
654
655 TProfile *p=0;
656 if (profx)
657 {
658 const TString pfx(Form("%sProfX", fHist->GetName()));
659
660 if (same && (p=dynamic_cast<TProfile*>(gPad->FindObject(pfx))))
661 *fLog << warn << "TProfile " << pfx << " already in pad." << endl;
662
663 p = ((TH2*)fHist)->ProfileX(pfx, -1, -1, "s");
664 p->UseCurrentStyle();
665 p->SetLineColor(blue ? kBlue : fHist->GetLineColor());
666 p->SetBit(kCanDelete);
667 p->SetDirectory(NULL);
668 p->SetXTitle(fHist->GetXaxis()->GetTitle());
669 p->SetYTitle(fHist->GetYaxis()->GetTitle());
670 p->Draw(only&&!same?"":"same");
671 }
672 if (profy)
673 {
674 const TString pfy(Form("%sProfY", fHist->GetName()));
675
676 if (same && (p=dynamic_cast<TProfile*>(gPad->FindObject(pfy))))
677 *fLog << warn << "TProfile " << pfy << " already in pad." << endl;
678
679 p = ((TH2*)fHist)->ProfileY(pfy, -1, -1, "s");
680 p->UseCurrentStyle();
681 p->SetLineColor(blue ? kBlue : fHist->GetLineColor());
682 p->SetBit(kCanDelete);
683 p->SetDirectory(NULL);
684 p->SetYTitle(fHist->GetXaxis()->GetTitle());
685 p->SetXTitle(fHist->GetYaxis()->GetTitle());
686 p->Draw(only&&!same?"":"same");
687 }
688
689 //AppendPad("log");
690}
691
692// --------------------------------------------------------------------------
693//
694// Implementation of SavePrimitive. Used to write the call to a constructor
695// to a macro. In the original root implementation it is used to write
696// gui elements to a macro-file.
697//
698void MH3::StreamPrimitive(ostream &out) const
699{
700 TString name = GetUniqueName();
701
702 out << " MH3 " << name << "(\"";
703 out << fData[0]->GetRule() << "\"";
704 if (fDimension>1)
705 out << ", \"" << fData[1]->GetRule() << "\"";
706 if (fDimension>2)
707 out << ", \"" << fData[2]->GetRule() << "\"";
708
709 out << ");" << endl;
710
711 if (fName!=gsDefName)
712 out << " " << name << ".SetName(\"" << fName << "\");" << endl;
713
714 if (fTitle!=gsDefTitle)
715 out << " " << name << ".SetTitle(\"" << fTitle << "\");" << endl;
716
717 switch (fDimension)
718 {
719 case 3:
720 if (fScale[2]!=1)
721 out << " " << name << ".SetScaleZ(" << fScale[2] << ");" << endl;
722 case 2:
723 if (fScale[1]!=1)
724 out << " " << name << ".SetScaleY(" << fScale[1] << ");" << endl;
725 case 1:
726 if (fScale[0]!=1)
727 out << " " << name << ".SetScaleX(" << fScale[0] << ");" << endl;
728 }
729}
730
731// --------------------------------------------------------------------------
732//
733// Used to rebuild a MH3 object of the same type (data members,
734// dimension, ...)
735//
736MParContainer *MH3::New() const
737{
738 MH3 *h = NULL;
739
740 if (fData[0] == NULL)
741 h=new MH3(fDimension);
742 else
743 switch (fDimension)
744 {
745 case 1:
746 h=new MH3(fData[0]->GetRule());
747 break;
748 case 2:
749 h=new MH3(fData[0]->GetRule(), fData[1]->GetRule());
750 break;
751 case 3:
752 h=new MH3(fData[0]->GetRule(), fData[1]->GetRule(), fData[2]->GetRule());
753 break;
754 }
755 switch (fDimension)
756 {
757 case 3:
758 h->SetScaleZ(fScale[2]);
759 case 2:
760 h->SetScaleY(fScale[1]);
761 case 1:
762 h->SetScaleX(fScale[0]);
763 }
764 return h;
765}
766
767TString MH3::GetRule(const Char_t axis) const
768{
769 switch (tolower(axis))
770 {
771 case 'x':
772 return fData[0] ? fData[0]->GetRule() : TString("");
773 case 'y':
774 return fData[1] ? fData[1]->GetRule() : TString("");
775 case 'z':
776 return fData[2] ? fData[2]->GetRule() : TString("");
777 default:
778 return "<n/a>";
779 }
780}
781
782// --------------------------------------------------------------------------
783//
784// Returns the total number of bins in a histogram (excluding under- and
785// overflow bins)
786//
787Int_t MH3::GetNbins() const
788{
789 Int_t num = 1;
790
791 switch (fDimension)
792 {
793 case 3:
794 num *= fHist->GetNbinsZ()+2;
795 case 2:
796 num *= fHist->GetNbinsY()+2;
797 case 1:
798 num *= fHist->GetNbinsX()+2;
799 }
800
801 return num;
802}
803
804// --------------------------------------------------------------------------
805//
806// Returns the total number of bins in a histogram (excluding under- and
807// overflow bins) Return -1 if bin is underflow or overflow
808//
809Int_t MH3::FindFixBin(Double_t x, Double_t y, Double_t z) const
810{
811 const TAxis &axex = *fHist->GetXaxis();
812 const TAxis &axey = *fHist->GetYaxis();
813 const TAxis &axez = *fHist->GetZaxis();
814
815 Int_t binz = 0;
816 Int_t biny = 0;
817 Int_t binx = 0;
818
819 switch (fDimension)
820 {
821 case 3:
822 binz = axez.FindFixBin(z);
823 if (binz>axez.GetLast() || binz<axez.GetFirst())
824 return -1;
825 case 2:
826 biny = axey.FindFixBin(y);
827 if (biny>axey.GetLast() || biny<axey.GetFirst())
828 return -1;
829 case 1:
830 binx = axex.FindFixBin(x);
831 if (binx<axex.GetFirst() || binx>axex.GetLast())
832 return -1;
833 }
834
835 const Int_t nx = fHist->GetNbinsX()+2;
836 const Int_t ny = fHist->GetNbinsY()+2;
837
838 return binx + nx*(biny +ny*binz);
839}
Note: See TracBrowser for help on using the repository browser.