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

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