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

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