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

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