Logo ROOT   6.10/09
Reference Guide
THnBase.cxx
Go to the documentation of this file.
1 // @(#)root/hist:$Id$
2 // Author: Axel Naumann (2011-12-20)
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2012, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 #include "THnBase.h"
13 
14 #include "TAxis.h"
15 #include "TBrowser.h"
16 #include "TBuffer.h"
17 #include "TError.h"
18 #include "TClass.h"
19 #include "TF1.h"
20 #include "TH1D.h"
21 #include "TH2D.h"
22 #include "TH3D.h"
23 #include "THn.h"
24 #include "THnSparse.h"
25 #include "TMath.h"
26 #include "TRandom.h"
27 #include "TVirtualPad.h"
28 
29 #include "HFitInterface.h"
30 #include "Fit/DataRange.h"
31 #include "Fit/SparseData.h"
32 #include "Math/MinimizerOptions.h"
33 #include "Math/WrappedMultiTF1.h"
34 
35 
36 /** \class THnBase
37  \ingroup Hist
38 Multidimensional histogram base.
39 Defines common functionality and interfaces for THn, THnSparse.
40 */
41 
43 
44 ////////////////////////////////////////////////////////////////////////////////
45 /// Construct a THnBase with "dim" dimensions,
46 /// "nbins" holds the number of bins for each dimension;
47 /// "xmin" and "xmax" the minimal and maximal value for each dimension.
48 /// The arrays "xmin" and "xmax" can be NULL; in that case SetBinEdges()
49 /// must be called for each dimension.
50 
51 THnBase::THnBase(const char* name, const char* title, Int_t dim,
52  const Int_t* nbins, const Double_t* xmin, const Double_t* xmax):
53 TNamed(name, title), fNdimensions(dim), fAxes(dim), fBrowsables(dim),
54 fEntries(0), fTsumw(0), fTsumw2(-1.), fTsumwx(dim), fTsumwx2(dim),
55 fIntegral(0), fIntegralStatus(kNoInt)
56 {
57  for (Int_t i = 0; i < fNdimensions; ++i) {
58  TAxis* axis = new TAxis(nbins[i], xmin ? xmin[i] : 0., xmax ? xmax[i] : 1.);
59  axis->SetName(TString::Format("axis%d", i));
60  fAxes.AddAtAndExpand(axis, i);
61  }
62  SetTitle(title);
63  fAxes.SetOwner();
64 }
65 
66 ////////////////////////////////////////////////////////////////////////////////
67 /// Destruct a THnBase
68 
70  if (fIntegralStatus != kNoInt) delete [] fIntegral;
71 }
72 
73 
74 ////////////////////////////////////////////////////////////////////////////////
75 /// Create a new THnBase object that is of the same type as *this,
76 /// but with dimensions and bins given by axes.
77 /// If keepTargetAxis is true, the axes will keep their original xmin / xmax,
78 /// else they will be restricted to the range selected (first / last).
79 
80 THnBase* THnBase::CloneEmpty(const char* name, const char* title,
81  const TObjArray* axes, Bool_t keepTargetAxis) const
82 {
83  THnBase* ret = (THnBase*)IsA()->New();
84  Int_t chunkSize = 1024 * 16;
86  chunkSize = ((const THnSparse*)this)->GetChunkSize();
87  }
88  ret->Init(name, title, axes, keepTargetAxis, chunkSize);
89  return ret;
90 }
91 
92 
93 ////////////////////////////////////////////////////////////////////////////////
94 /// Initialize axes and name.
95 
96 void THnBase::Init(const char* name, const char* title,
97  const TObjArray* axes, Bool_t keepTargetAxis,
98  Int_t chunkSize /*= 1024 * 16*/)
99 {
100  SetNameTitle(name, title);
101 
102  TIter iAxis(axes);
103  const TAxis* axis = 0;
104  Int_t pos = 0;
105  Int_t *nbins = new Int_t[axes->GetEntriesFast()];
106  while ((axis = (TAxis*)iAxis())) {
107  TAxis* reqaxis = new TAxis(*axis);
108  if (!keepTargetAxis && axis->TestBit(TAxis::kAxisRange)) {
109  Int_t binFirst = axis->GetFirst();
110  // The lowest egde of the underflow is meaningless.
111  if (binFirst == 0)
112  binFirst = 1;
113  Int_t binLast = axis->GetLast();
114  // The overflow edge is implicit.
115  if (binLast > axis->GetNbins())
116  binLast = axis->GetNbins();
117  Int_t nBins = binLast - binFirst + 1;
118  if (axis->GetXbins()->GetSize()) {
119  // non-uniform bins:
120  reqaxis->Set(nBins, axis->GetXbins()->GetArray() + binFirst - 1);
121  } else {
122  // uniform bins:
123  reqaxis->Set(nBins, axis->GetBinLowEdge(binFirst), axis->GetBinUpEdge(binLast));
124  }
125  reqaxis->ResetBit(TAxis::kAxisRange);
126  }
127 
128  nbins[pos] = reqaxis->GetNbins();
129  fAxes.AddAtAndExpand(new TAxis(*reqaxis), pos++);
130  }
131  fAxes.SetOwner();
132 
133  fNdimensions = axes->GetEntriesFast();
134  InitStorage(nbins, chunkSize);
135  delete [] nbins;
136 }
137 
138 
139 ////////////////////////////////////////////////////////////////////////////////
140 /// Create an empty histogram with name and title with a given
141 /// set of axes. Create a TH1D/TH2D/TH3D, depending on the number
142 /// of elements in axes.
143 
144 TH1* THnBase::CreateHist(const char* name, const char* title,
145  const TObjArray* axes,
146  Bool_t keepTargetAxis ) const {
147  const int ndim = axes->GetSize();
148 
149  TH1* hist = 0;
150  // create hist with dummy axes, we fix them later.
151  if (ndim == 1)
152  hist = new TH1D(name, title, 1, 0., 1.);
153  else if (ndim == 2)
154  hist = new TH2D(name, title, 1, 0., 1., 1, 0., 1.);
155  else if (ndim == 3)
156  hist = new TH3D(name, title, 1, 0., 1., 1, 0., 1., 1, 0., 1.);
157  else {
158  Error("CreateHist", "Cannot create histogram %s with %d dimensions!", name, ndim);
159  return 0;
160  }
161 
162  TAxis* hax[3] = {hist->GetXaxis(), hist->GetYaxis(), hist->GetZaxis()};
163  for (Int_t d = 0; d < ndim; ++d) {
164  TAxis* reqaxis = (TAxis*)(*axes)[d];
165  hax[d]->SetTitle(reqaxis->GetTitle());
166  if (!keepTargetAxis && reqaxis->TestBit(TAxis::kAxisRange)) {
167  // axis cannot extend to underflow/overlflows (fix ROOT-8781)
168  Int_t binFirst = std::max(reqaxis->GetFirst(),1);
169  Int_t binLast = std::min(reqaxis->GetLast(), reqaxis->GetNbins() );
170  Int_t nBins = binLast - binFirst + 1;
171  if (reqaxis->GetXbins()->GetSize()) {
172  // non-uniform bins:
173  hax[d]->Set(nBins, reqaxis->GetXbins()->GetArray() + binFirst - 1);
174  } else {
175  // uniform bins:
176  hax[d]->Set(nBins, reqaxis->GetBinLowEdge(binFirst), reqaxis->GetBinUpEdge(binLast));
177  }
178  } else {
179  if (reqaxis->GetXbins()->GetSize()) {
180  // non-uniform bins:
181  hax[d]->Set(reqaxis->GetNbins(), reqaxis->GetXbins()->GetArray());
182  } else {
183  // uniform bins:
184  hax[d]->Set(reqaxis->GetNbins(), reqaxis->GetXmin(), reqaxis->GetXmax());
185  }
186  }
187  }
188 
189  hist->Rebuild();
190 
191  return hist;
192 }
193 
194 ////////////////////////////////////////////////////////////////////////////////
195 /// Create a THn / THnSparse object from a histogram deriving from TH1.
196 
197 THnBase* THnBase::CreateHnAny(const char* name, const char* title,
198  const TH1* h, Bool_t sparse, Int_t chunkSize)
199 {
200  // Get the dimension of the TH1
201  int ndim = h->GetDimension();
202 
203  // Axis properties
204  int nbins[3] = {0,0,0};
205  double minRange[3] = {0.,0.,0.};
206  double maxRange[3] = {0.,0.,0.};
207  const TAxis* axis[3] = { h->GetXaxis(), h->GetYaxis(), h->GetZaxis() };
208  for (int i = 0; i < ndim; ++i) {
209  nbins[i] = axis[i]->GetNbins();
210  minRange[i] = axis[i]->GetXmin();
211  maxRange[i] = axis[i]->GetXmax();
212  }
213 
214  // Create the corresponding THnSparse, depending on the storage
215  // type of the TH1. The class name will be "TH??\0" where the first
216  // ? is 1,2 or 3 and the second ? indicates the stograge as C, S,
217  // I, F or D.
218  THnBase* s = 0;
219  const char* cname( h->ClassName() );
220  if (cname[0] == 'T' && cname[1] == 'H'
221  && cname[2] >= '1' && cname[2] <= '3' && cname[4] == 0) {
222 
223 #define R__THNBCASE(TAG) \
224 if (sparse) { \
225 s = new _NAME2_(THnSparse,TAG)(name, title, ndim, nbins, \
226 minRange, maxRange, chunkSize); \
227 } else { \
228 s = new _NAME2_(THn,TAG)(name, title, ndim, nbins, \
229 minRange, maxRange); \
230 } \
231 break;
232 
233  switch (cname[3]) {
234  case 'F': R__THNBCASE(F);
235  case 'D': R__THNBCASE(D);
236  case 'I': R__THNBCASE(I);
237  case 'S': R__THNBCASE(S);
238  case 'C': R__THNBCASE(C);
239  }
240 #undef R__THNBCASE
241  }
242  if (!s) {
243  ::Warning("THnSparse::CreateHnAny", "Unknown Type of Histogram");
244  return 0;
245  }
246 
247  for (int i = 0; i < ndim; ++i) {
248  s->GetAxis(i)->SetTitle(axis[i]->GetTitle());
249  }
250 
251  // Get the array to know the number of entries of the TH1
252  const TArray *array = dynamic_cast<const TArray*>(h);
253  if (!array) {
254  ::Warning("THnSparse::CreateHnAny", "Unknown Type of Histogram");
255  return 0;
256  }
257 
258  s->Add(h);
259  return s;
260 }
261 
262 
263 ////////////////////////////////////////////////////////////////////////////////
264 /// Create a THnSparse (if "sparse") or THn from "hn", possibly
265 /// converting THn <-> THnSparse.
266 
267 THnBase* THnBase::CreateHnAny(const char* name, const char* title,
268  const THnBase* hn, Bool_t sparse,
269  Int_t chunkSize /*= 1024 * 16*/)
270 {
271  TClass* type = 0;
272  if (hn->InheritsFrom(THnSparse::Class())) {
273  if (sparse) type = hn->IsA();
274  else {
275  char bintype;
276  if (hn->InheritsFrom(THnSparseD::Class())) bintype = 'D';
277  else if (hn->InheritsFrom(THnSparseF::Class())) bintype = 'F';
278  else if (hn->InheritsFrom(THnSparseL::Class())) bintype = 'L';
279  else if (hn->InheritsFrom(THnSparseI::Class())) bintype = 'I';
280  else if (hn->InheritsFrom(THnSparseS::Class())) bintype = 'S';
281  else if (hn->InheritsFrom(THnSparseC::Class())) bintype = 'C';
282  else {
283  hn->Error("CreateHnAny", "Type %s not implemented; please inform the ROOT team!",
284  hn->IsA()->GetName());
285  return 0;
286  }
287  type = TClass::GetClass(TString::Format("THn%c", bintype));
288  }
289  } else if (hn->InheritsFrom(THn::Class())) {
290  if (!sparse) type = hn->IsA();
291  else {
292  char bintype = 0;
293  if (hn->InheritsFrom(THnD::Class())) bintype = 'D';
294  else if (hn->InheritsFrom(THnF::Class())) bintype = 'F';
295  else if (hn->InheritsFrom(THnC::Class())) bintype = 'C';
296  else if (hn->InheritsFrom(THnS::Class())) bintype = 'S';
297  else if (hn->InheritsFrom(THnI::Class())) bintype = 'I';
298  else if (hn->InheritsFrom(THnL::Class())) bintype = 'L';
299  else if (hn->InheritsFrom(THnL64::Class())) {
300  hn->Error("CreateHnAny", "Type THnSparse with Long64_t bins is not available!");
301  return 0;
302  }
303  if (bintype) {
304  type = TClass::GetClass(TString::Format("THnSparse%c", bintype));
305  }
306  }
307  } else {
308  hn->Error("CreateHnAny", "Unhandled type %s, not deriving from THn nor THnSparse!",
309  hn->IsA()->GetName());
310  return 0;
311  }
312  if (!type) {
313  hn->Error("CreateHnAny", "Unhandled type %s, please inform the ROOT team!",
314  hn->IsA()->GetName());
315  return 0;
316  }
317 
318  THnBase* ret = (THnBase*)type->New();
319  ret->Init(name, title, hn->GetListOfAxes(),
320  kFALSE /*keepTargetAxes*/, chunkSize);
321 
322  ret->Add(hn);
323  return ret;
324 }
325 
326 
327 ////////////////////////////////////////////////////////////////////////////////
328 /// Fill the THnBase with the bins of hist that have content
329 /// or error != 0.
330 
331 void THnBase::Add(const TH1* hist, Double_t c /*=1.*/)
332 {
333  Long64_t nbins = hist->GetNcells();
334  int x[3] = {0,0,0};
335  for (int i = 0; i < nbins; ++i) {
336  double value = hist->GetBinContent(i);
337  double error = hist->GetBinError(i);
338  if (!value && !error) continue;
339  hist->GetBinXYZ(i, x[0], x[1], x[2]);
340  SetBinContent(x, value * c);
341  SetBinError(x, error * c);
342  }
343 }
344 
345 ////////////////////////////////////////////////////////////////////////////////
346 /// Fit a THnSparse with function f
347 ///
348 /// since the data is sparse by default a likelihood fit is performed
349 /// merging all the regions with empty bins for betetr performance efficiency
350 ///
351 /// Since the THnSparse is not drawn no graphics options are passed
352 /// Here is the list of possible options
353 ///
354 /// = "I" Use integral of function in bin instead of value at bin center
355 /// = "X" Use chi2 method (default is log-likelihood method)
356 /// = "U" Use a User specified fitting algorithm (via SetFCN)
357 /// = "Q" Quiet mode (minimum printing)
358 /// = "V" Verbose mode (default is between Q and V)
359 /// = "E" Perform better Errors estimation using Minos technique
360 /// = "B" Use this option when you want to fix one or more parameters
361 /// and the fitting function is like "gaus", "expo", "poln", "landau".
362 /// = "M" More. Improve fit results
363 /// = "R" Use the Range specified in the function range
364 
366 {
367 
368  Foption_t fitOption;
369 
370  if (!TH1::FitOptionsMake(option,fitOption)) return 0;
371 
372  // The function used to fit cannot be stored in a THnSparse. It
373  // cannot be drawn either. Perhaps in the future.
374  fitOption.Nostore = true;
375  // Use likelihood fit if not specified
376  if (!fitOption.Chi2) fitOption.Like = true;
377  // create range and minimizer options with default values
379  for ( int i = 0; i < GetNdimensions(); ++i ) {
380  TAxis *axis = GetAxis(i);
381  range.AddRange(i, axis->GetXmin(), axis->GetXmax());
382  }
384 
385  return ROOT::Fit::FitObject(this, f , fitOption , minOption, goption, range);
386 }
387 
388 ////////////////////////////////////////////////////////////////////////////////
389 /// Generate an n-dimensional random tuple based on the histogrammed
390 /// distribution. If subBinRandom, the returned tuple will be additionally
391 /// randomly distributed within the randomized bin, using a flat
392 /// distribution.
393 
394 void THnBase::GetRandom(Double_t *rand, Bool_t subBinRandom /* = kTRUE */)
395 {
396  // check whether the integral array is valid
397  if (fIntegralStatus != kValidInt)
398  ComputeIntegral();
399 
400  // generate a random bin
401  Double_t p = gRandom->Rndm();
403  const Int_t nStaticBins = 40;
404  Int_t bin[nStaticBins];
405  Int_t* pBin = bin;
406  if (GetNdimensions() > nStaticBins) {
407  pBin = new Int_t[GetNdimensions()];
408  }
409  GetBinContent(idx, pBin);
410 
411  // convert bin coordinates to real values
412  for (Int_t i = 0; i < fNdimensions; i++) {
413  rand[i] = GetAxis(i)->GetBinCenter(pBin[i]);
414 
415  // randomize the vector withing a bin
416  if (subBinRandom)
417  rand[i] += (gRandom->Rndm() - 0.5) * GetAxis(i)->GetBinWidth(pBin[i]);
418  }
419  if (pBin != bin) {
420  delete [] pBin;
421  }
422 
423  return;
424 }
425 
426 ////////////////////////////////////////////////////////////////////////////////
427 /// Check whether bin coord is in range, as defined by TAxis::SetRange().
428 
430 {
431  Int_t min = 0;
432  Int_t max = 0;
433  for (Int_t i = 0; i < fNdimensions; ++i) {
434  TAxis *axis = GetAxis(i);
435  if (!axis->TestBit(TAxis::kAxisRange)) continue;
436  min = axis->GetFirst();
437  max = axis->GetLast();
438  if (coord[i] < min || coord[i] > max)
439  return kFALSE;
440  }
441  return kTRUE;
442 }
443 
444 ////////////////////////////////////////////////////////////////////////////////
445 /// Project all bins into a ndim-dimensional THn / THnSparse (whatever
446 /// *this is) or if (ndim < 4 and !wantNDim) a TH1/2/3 histogram,
447 /// keeping only axes in dim (specifying ndim dimensions).
448 /// If "option" contains "E" errors will be calculated.
449 /// "A" ranges of the target axes will be ignored.
450 /// "O" original axis range of the target axes will be
451 /// kept, but only bins inside the selected range
452 /// will be filled.
453 
455  Bool_t wantNDim,
456  Option_t* option /*= ""*/) const
457 {
458  TString name(GetName());
459  name +="_proj";
460 
461  for (Int_t d = 0; d < ndim; ++d) {
462  name += "_";
463  name += dim[d];
464  }
465 
466  TString title(GetTitle());
467  Ssiz_t posInsert = title.First(';');
468  if (posInsert == kNPOS) {
469  title += " projection ";
470  for (Int_t d = 0; d < ndim; ++d)
471  title += GetAxis(dim[d])->GetTitle();
472  } else {
473  for (Int_t d = ndim - 1; d >= 0; --d) {
474  title.Insert(posInsert, GetAxis(d)->GetTitle());
475  if (d)
476  title.Insert(posInsert, ", ");
477  }
478  title.Insert(posInsert, " projection ");
479  }
480 
481  TObjArray newaxes(ndim);
482  for (Int_t d = 0; d < ndim; ++d) {
483  newaxes.AddAt(GetAxis(dim[d]),d);
484  }
485 
486  THnBase* hn = 0;
487  TH1* hist = 0;
488  TObject* ret = 0;
489 
490  Bool_t* hadRange = 0;
491  Bool_t ignoreTargetRange = (option && (strchr(option, 'A') || strchr(option, 'a')));
492  Bool_t keepTargetAxis = ignoreTargetRange || (option && (strchr(option, 'O') || strchr(option, 'o')));
493  if (ignoreTargetRange) {
494  hadRange = new Bool_t[ndim];
495  for (Int_t d = 0; d < ndim; ++d){
496  TAxis *axis = GetAxis(dim[d]);
497  hadRange[d] = axis->TestBit(TAxis::kAxisRange);
499  }
500  }
501 
502  if (wantNDim)
503  ret = hn = CloneEmpty(name, title, &newaxes, keepTargetAxis);
504  else
505  ret = hist = CreateHist(name, title, &newaxes, keepTargetAxis);
506 
507  if (keepTargetAxis) {
508  // make the whole axes visible, i.e. unset the range
509  if (wantNDim) {
510  for (Int_t d = 0; d < ndim; ++d) {
511  hn->GetAxis(d)->SetRange(0, 0);
512  }
513  } else {
514  hist->GetXaxis()->SetRange(0, 0);
515  hist->GetYaxis()->SetRange(0, 0);
516  hist->GetZaxis()->SetRange(0, 0);
517  }
518  }
519 
520  Bool_t haveErrors = GetCalculateErrors();
521  Bool_t wantErrors = haveErrors || (option && (strchr(option, 'E') || strchr(option, 'e')));
522 
523  Int_t* bins = new Int_t[ndim];
524  Long64_t myLinBin = 0;
525 
526  THnIter iter(this, kTRUE /*use axis range*/);
527 
528  while ((myLinBin = iter.Next()) >= 0) {
529  Double_t v = GetBinContent(myLinBin);
530 
531  for (Int_t d = 0; d < ndim; ++d) {
532  bins[d] = iter.GetCoord(dim[d]);
533  if (!keepTargetAxis && GetAxis(dim[d])->TestBit(TAxis::kAxisRange)) {
534  Int_t binOffset = GetAxis(dim[d])->GetFirst();
535  // Don't subtract even more if underflow is alreday included:
536  if (binOffset > 0) --binOffset;
537  bins[d] -= binOffset;
538  }
539  }
540 
541  Long64_t targetLinBin = -1;
542  if (!wantNDim) {
543  if (ndim == 1) targetLinBin = bins[0];
544  else if (ndim == 2) targetLinBin = hist->GetBin(bins[0], bins[1]);
545  else if (ndim == 3) targetLinBin = hist->GetBin(bins[0], bins[1], bins[2]);
546  } else {
547  targetLinBin = hn->GetBin(bins, kTRUE /*allocate*/);
548  }
549 
550  if (wantErrors) {
551  Double_t err2 = 0.;
552  if (haveErrors) {
553  err2 = GetBinError2(myLinBin);
554  } else {
555  err2 = v;
556  }
557  if (wantNDim) {
558  hn->AddBinError2(targetLinBin, err2);
559  } else {
560  Double_t preverr = hist->GetBinError(targetLinBin);
561  hist->SetBinError(targetLinBin, TMath::Sqrt(preverr * preverr + err2));
562  }
563  }
564 
565  // only _after_ error calculation, or sqrt(v) is taken into account!
566  if (wantNDim)
567  hn->AddBinContent(targetLinBin, v);
568  else
569  hist->AddBinContent(targetLinBin, v);
570  }
571 
572  delete [] bins;
573 
574  if (wantNDim) {
575  hn->SetEntries(fEntries);
576  } else {
577  if (!iter.HaveSkippedBin()) {
578  hist->SetEntries(fEntries);
579  } else {
580  // re-compute the entries
581  // in case of error calculation (i.e. when Sumw2() is set)
582  // use the effective entries for the entries
583  // since this is the only way to estimate them
584  hist->ResetStats();
585  Double_t entries = hist->GetEffectiveEntries();
586  if (!wantErrors) {
587  // to avoid numerical rounding
588  entries = TMath::Floor(entries + 0.5);
589  }
590  hist->SetEntries(entries);
591  }
592  }
593 
594  if (hadRange) {
595  // reset kAxisRange bit:
596  for (Int_t d = 0; d < ndim; ++d)
597  GetAxis(dim[d])->SetBit(TAxis::kAxisRange, hadRange[d]);
598 
599  delete [] hadRange;
600  }
601 
602  return ret;
603 }
604 
605 ////////////////////////////////////////////////////////////////////////////////
606 /// Scale contents and errors of this histogram by c:
607 /// this = this * c
608 /// It does not modify the histogram's number of entries.
609 
611 {
612 
613  Double_t nEntries = GetEntries();
614  // Scale the contents & errors
615  Bool_t haveErrors = GetCalculateErrors();
616  Long64_t i = 0;
617  THnIter iter(this);
618  while ((i = iter.Next()) >= 0) {
619  // Get the content of the bin from the current histogram
620  Double_t v = GetBinContent(i);
621  SetBinContent(i, c * v);
622  if (haveErrors) {
623  Double_t err2 = GetBinError2(i);
624  SetBinError2(i, c * c * err2);
625  }
626  }
627  SetEntries(nEntries);
628 }
629 
630 ////////////////////////////////////////////////////////////////////////////////
631 /// Add() implementation for both rebinned histograms and those with identical
632 /// binning. See THnBase::Add().
633 
634 void THnBase::AddInternal(const THnBase* h, Double_t c, Bool_t rebinned)
635 {
636  if (fNdimensions != h->GetNdimensions()) {
637  Warning("RebinnedAdd", "Different number of dimensions, cannot carry out operation on the histograms");
638  return;
639  }
640 
641  // Trigger error calculation if h has it
642  if (!GetCalculateErrors() && h->GetCalculateErrors())
643  Sumw2();
644  Bool_t haveErrors = GetCalculateErrors();
645 
646  Double_t* x = 0;
647  if (rebinned) {
648  x = new Double_t[fNdimensions];
649  }
650  Int_t* coord = new Int_t[fNdimensions];
651 
652  // Expand the exmap if needed, to reduce collisions
653  Long64_t numTargetBins = GetNbins() + h->GetNbins();
654  Reserve(numTargetBins);
655 
656  Long64_t i = 0;
657  THnIter iter(h);
658  // Add to this whatever is found inside the other histogram
659  while ((i = iter.Next(coord)) >= 0) {
660  // Get the content of the bin from the second histogram
661  Double_t v = h->GetBinContent(i);
662 
663  Long64_t mybinidx = -1;
664  if (rebinned) {
665  // Get the bin center given a coord
666  for (Int_t j = 0; j < fNdimensions; ++j)
667  x[j] = h->GetAxis(j)->GetBinCenter(coord[j]);
668 
669  mybinidx = GetBin(x, kTRUE /* allocate*/);
670  } else {
671  mybinidx = GetBin(coord, kTRUE /*allocate*/);
672  }
673 
674  if (haveErrors) {
675  Double_t err2 = h->GetBinError2(i) * c * c;
676  AddBinError2(mybinidx, err2);
677  }
678  // only _after_ error calculation, or sqrt(v) is taken into account!
679  AddBinContent(mybinidx, c * v);
680  }
681 
682  delete [] coord;
683  delete [] x;
684 
685  Double_t nEntries = GetEntries() + c * h->GetEntries();
686  SetEntries(nEntries);
687 }
688 
689 ////////////////////////////////////////////////////////////////////////////////
690 /// Add contents of h scaled by c to this histogram:
691 /// this = this + c * h
692 /// Note that if h has Sumw2 set, Sumw2 is automatically called for this
693 /// if not already set.
694 
696 {
697  // Check consistency of the input
698  if (!CheckConsistency(h, "Add")) return;
699 
700  AddInternal(h, c, kFALSE);
701 }
702 
703 ////////////////////////////////////////////////////////////////////////////////
704 /// Add contents of h scaled by c to this histogram:
705 /// this = this + c * h
706 /// Note that if h has Sumw2 set, Sumw2 is automatically called for this
707 /// if not already set.
708 /// In contrast to Add(), RebinnedAdd() does not require consistent binning of
709 /// this and h; instead, each bin's center is used to determine the target bin.
710 
712 {
713  AddInternal(h, c, kTRUE);
714 }
715 
716 
717 ////////////////////////////////////////////////////////////////////////////////
718 /// Merge this with a list of THnBase's. All THnBase's provided
719 /// in the list must have the same bin layout!
720 
722 {
723  if (!list) return 0;
724  if (list->IsEmpty()) return (Long64_t)GetEntries();
725 
726  Long64_t sumNbins = GetNbins();
727  TIter iter(list);
728  const TObject* addMeObj = 0;
729  while ((addMeObj = iter())) {
730  const THnBase* addMe = dynamic_cast<const THnBase*>(addMeObj);
731  if (addMe) {
732  sumNbins += addMe->GetNbins();
733  }
734  }
735  Reserve(sumNbins);
736 
737  iter.Reset();
738  while ((addMeObj = iter())) {
739  const THnBase* addMe = dynamic_cast<const THnBase*>(addMeObj);
740  if (!addMe)
741  Error("Merge", "Object named %s is not THnBase! Skipping it.",
742  addMeObj->GetName());
743  else
744  Add(addMe);
745  }
746  return (Long64_t)GetEntries();
747 }
748 
749 
750 ////////////////////////////////////////////////////////////////////////////////
751 /// Multiply this histogram by histogram h
752 /// this = this * h
753 /// Note that if h has Sumw2 set, Sumw2 is automatically called for this
754 /// if not already set.
755 
757 {
758  // Check consistency of the input
759  if(!CheckConsistency(h, "Multiply"))return;
760 
761  // Trigger error calculation if h has it
762  Bool_t wantErrors = kFALSE;
764  wantErrors = kTRUE;
765 
766  if (wantErrors) Sumw2();
767 
768  Double_t nEntries = GetEntries();
769  // Now multiply the contents: in this case we have the intersection of the sets of bins
770  Int_t* coord = new Int_t[fNdimensions];
771  Long64_t i = 0;
772  THnIter iter(this);
773  // Add to this whatever is found inside the other histogram
774  while ((i = iter.Next(coord)) >= 0) {
775  // Get the content of the bin from the current histogram
776  Double_t v1 = GetBinContent(i);
777  // Now look at the bin with the same coordinates in h
778  Long64_t idxh = h->GetBin(coord);
779  Double_t v2 = 0.;
780  if (idxh >= 0) v2 = h->GetBinContent(idxh);
781  SetBinContent(i, v1 * v2);
782  if (wantErrors) {
783  Double_t err1 = GetBinError(i) * v2;
784  Double_t err2 = 0.;
785  if (idxh >= 0) err2 = h->GetBinError(idxh) * v1;
786  SetBinError(i, TMath::Sqrt((err2 * err2 + err1 * err1)));
787  }
788  }
789  SetEntries(nEntries);
790 
791  delete [] coord;
792 }
793 
794 ////////////////////////////////////////////////////////////////////////////////
795 /// Performs the operation: this = this*c*f1
796 /// if errors are defined, errors are also recalculated.
797 ///
798 /// Only bins inside the function range are recomputed.
799 /// IMPORTANT NOTE: If you intend to use the errors of this histogram later
800 /// you should call Sumw2 before making this operation.
801 /// This is particularly important if you fit the histogram after
802 /// calling Multiply()
803 
805 {
806  Int_t* coord = new Int_t[fNdimensions];
808 
809  Bool_t wantErrors = GetCalculateErrors();
810  if (wantErrors) Sumw2();
811 
812  Long64_t i = 0;
813  THnIter iter(this);
814  // Add to this whatever is found inside the other histogram
815  while ((i = iter.Next(coord)) >= 0) {
816  Double_t value = GetBinContent(i);
817 
818  // Get the bin coordinates given an index array
819  for (Int_t j = 0; j < fNdimensions; ++j)
820  x[j] = GetAxis(j)->GetBinCenter(coord[j]);
821 
822  if (!f->IsInside(x))
823  continue;
825 
826  // Evaluate function at points
827  Double_t fvalue = f->EvalPar(x, NULL);
828 
829  SetBinContent(i, c * fvalue * value);
830  if (wantErrors) {
831  Double_t error = GetBinError(i);
832  SetBinError(i, c * fvalue * error);
833  }
834  }
835 
836  delete [] x;
837  delete [] coord;
838 }
839 
840 ////////////////////////////////////////////////////////////////////////////////
841 /// Divide this histogram by h
842 /// this = this/(h)
843 /// Note that if h has Sumw2 set, Sumw2 is automatically called for
844 /// this if not already set.
845 /// The resulting errors are calculated assuming uncorrelated content.
846 
848 {
849  // Check consistency of the input
850  if (!CheckConsistency(h, "Divide"))return;
851 
852  // Trigger error calculation if h has it
853  Bool_t wantErrors=GetCalculateErrors();
854  if (!GetCalculateErrors() && h->GetCalculateErrors())
855  wantErrors=kTRUE;
856 
857  // Remember original histogram statistics
858  Double_t nEntries = fEntries;
859 
860  if (wantErrors) Sumw2();
861  Bool_t didWarn = kFALSE;
862 
863  // Now divide the contents: also in this case we have the intersection of the sets of bins
864  Int_t* coord = new Int_t[fNdimensions];
865  Long64_t i = 0;
866  THnIter iter(this);
867  // Add to this whatever is found inside the other histogram
868  while ((i = iter.Next(coord)) >= 0) {
869  // Get the content of the bin from the first histogram
870  Double_t v1 = GetBinContent(i);
871  // Now look at the bin with the same coordinates in h
872  Long64_t hbin = h->GetBin(coord);
873  Double_t v2 = h->GetBinContent(hbin);
874  if (!v2) {
875  v1 = 0.;
876  v2 = 1.;
877  if (!didWarn) {
878  Warning("Divide(h)", "Histogram h has empty bins - division by zero! Setting bin to 0.");
879  didWarn = kTRUE;
880  }
881  }
882  SetBinContent(i, v1 / v2);
883  if (wantErrors) {
884  Double_t err1 = GetBinError(i) * v2;
885  Double_t err2 = h->GetBinError(hbin) * v1;
886  Double_t b22 = v2 * v2;
887  Double_t err = (err1 * err1 + err2 * err2) / (b22 * b22);
888  SetBinError2(i, err);
889  }
890  }
891  delete [] coord;
892  SetEntries(nEntries);
893 }
894 
895 ////////////////////////////////////////////////////////////////////////////////
896 /// Replace contents of this histogram by multiplication of h1 by h2
897 /// this = (c1*h1)/(c2*h2)
898 /// Note that if h1 or h2 have Sumw2 set, Sumw2 is automatically called for
899 /// this if not already set.
900 /// The resulting errors are calculated assuming uncorrelated content.
901 /// However, if option ="B" is specified, Binomial errors are computed.
902 /// In this case c1 and c2 do not make real sense and they are ignored.
903 
904 void THnBase::Divide(const THnBase *h1, const THnBase *h2, Double_t c1, Double_t c2, Option_t *option)
905 {
906 
907  TString opt = option;
908  opt.ToLower();
909  Bool_t binomial = kFALSE;
910  if (opt.Contains("b")) binomial = kTRUE;
911 
912  // Check consistency of the input
913  if (!CheckConsistency(h1, "Divide") || !CheckConsistency(h2, "Divide"))return;
914  if (!c2) {
915  Error("Divide","Coefficient of dividing histogram cannot be zero");
916  return;
917  }
918 
919  Reset();
920 
921  // Trigger error calculation if h1 or h2 have it
922  if (!GetCalculateErrors() && (h1->GetCalculateErrors()|| h2->GetCalculateErrors() != 0))
923  Sumw2();
924 
925  // Count filled bins
926  Long64_t nFilledBins=0;
927 
928  // Now divide the contents: we have the intersection of the sets of bins
929 
930  Int_t* coord = new Int_t[fNdimensions];
931  memset(coord, 0, sizeof(Int_t) * fNdimensions);
932  Bool_t didWarn = kFALSE;
933 
934  Long64_t i = 0;
935  THnIter iter(h1);
936  // Add to this whatever is found inside the other histogram
937  while ((i = iter.Next(coord)) >= 0) {
938  // Get the content of the bin from the first histogram
939  Double_t v1 = h1->GetBinContent(i);
940  // Now look at the bin with the same coordinates in h2
941  Long64_t h2bin = h2->GetBin(coord);
942  Double_t v2 = h2->GetBinContent(h2bin);
943  if (!v2) {
944  v1 = 0.;
945  v2 = 1.;
946  if (!didWarn) {
947  Warning("Divide(h1, h2)", "Histogram h2 has empty bins - division by zero! Setting bin to 0.");
948  didWarn = kTRUE;
949  }
950  }
951  nFilledBins++;
952  Long64_t myBin = GetBin(coord);
953  SetBinContent(myBin, c1 * v1 / c2 / v2);
954  if(GetCalculateErrors()){
955  Double_t err1 = h1->GetBinError(i);
956  Double_t err2 = h2->GetBinError(h2bin);
957  Double_t errSq = 0.;
958  if (binomial) {
959  if (v1 != v2) {
960  Double_t w = v1 / v2;
961  err2 *= w;
962  errSq = TMath::Abs( ( (1. - 2.*w) * err1 * err1 + err2 * err2 ) / (v2 * v2) );
963  }
964  } else {
965  c1 *= c1;
966  c2 *= c2;
967  Double_t b22 = v2 * v2 * c2;
968  err1 *= v2;
969  err2 *= v1;
970  errSq = c1 * c2 * (err1 * err1 + err2 * err2) / (b22 * b22);
971  }
972  SetBinError2(myBin, errSq);
973  }
974  }
975 
976  delete [] coord;
977  SetFilledBins(nFilledBins);
978 
979  // Set as entries in the result histogram the entries in the numerator
980  SetEntries(h1->GetEntries());
981 }
982 
983 ////////////////////////////////////////////////////////////////////////////////
984 /// Consistency check on (some of) the parameters of two histograms (for operations).
985 
986 Bool_t THnBase::CheckConsistency(const THnBase *h, const char *tag) const
987 {
988  if (fNdimensions != h->GetNdimensions()) {
989  Warning(tag, "Different number of dimensions, cannot carry out operation on the histograms");
990  return kFALSE;
991  }
992  for (Int_t dim = 0; dim < fNdimensions; dim++){
993  if (GetAxis(dim)->GetNbins() != h->GetAxis(dim)->GetNbins()) {
994  Warning(tag, "Different number of bins on axis %i, cannot carry out operation on the histograms", dim);
995  return kFALSE;
996  }
997  }
998  return kTRUE;
999 }
1000 
1001 ////////////////////////////////////////////////////////////////////////////////
1002 /// Set the axis # of bins and bin limits on dimension idim
1003 
1004 void THnBase::SetBinEdges(Int_t idim, const Double_t* bins)
1005 {
1006  TAxis* axis = (TAxis*) fAxes[idim];
1007  axis->Set(axis->GetNbins(), bins);
1008 }
1009 
1010 ////////////////////////////////////////////////////////////////////////////////
1011 /// Change (i.e. set) the title.
1012 ///
1013 /// If title is in the form "stringt;string0;string1;string2 ..."
1014 /// the histogram title is set to stringt, the title of axis0 to string0,
1015 /// of axis1 to string1, of axis2 to string2, etc, just like it is done
1016 /// for TH1/TH2/TH3.
1017 /// To insert the character ";" in one of the titles, one should use "#;"
1018 /// or "#semicolon".
1019 
1020 void THnBase::SetTitle(const char *title)
1021 {
1022  fTitle = title;
1023  fTitle.ReplaceAll("#;",2,"#semicolon",10);
1024 
1025  Int_t endHistTitle = fTitle.First(';');
1026  if (endHistTitle >= 0) {
1027  // title contains a ';' so parse the axis titles
1028  Int_t posTitle = endHistTitle + 1;
1029  Int_t lenTitle = fTitle.Length();
1030  Int_t dim = 0;
1031  while (posTitle > 0 && posTitle < lenTitle && dim < fNdimensions){
1032  Int_t endTitle = fTitle.Index(";", posTitle);
1033  TString axisTitle = fTitle(posTitle, endTitle - posTitle);
1034  axisTitle.ReplaceAll("#semicolon", 10, ";", 1);
1035  GetAxis(dim)->SetTitle(axisTitle);
1036  dim++;
1037  if (endTitle > 0)
1038  posTitle = endTitle + 1;
1039  else
1040  posTitle = -1;
1041  }
1042  // Remove axis titles from histogram title
1043  fTitle.Remove(endHistTitle, lenTitle - endHistTitle);
1044  }
1045 
1046  fTitle.ReplaceAll("#semicolon", 10, ";", 1);
1047 
1048 }
1049 
1050 ////////////////////////////////////////////////////////////////////////////////
1051 /// Combine the content of "group" neighboring bins into
1052 /// a new bin and return the resulting THnBase.
1053 /// For group=2 and a 3 dimensional histogram, all "blocks"
1054 /// of 2*2*2 bins will be put into a bin.
1055 
1057 {
1058  Int_t* ngroup = new Int_t[GetNdimensions()];
1059  for (Int_t d = 0; d < GetNdimensions(); ++d)
1060  ngroup[d] = group;
1061  THnBase* ret = RebinBase(ngroup);
1062  delete [] ngroup;
1063  return ret;
1064 }
1065 
1066 ////////////////////////////////////////////////////////////////////////////////
1067 /// Combine the content of "group" neighboring bins for each dimension
1068 /// into a new bin and return the resulting THnBase.
1069 /// For group={2,1,1} and a 3 dimensional histogram, pairs of x-bins
1070 /// will be grouped.
1071 
1073 {
1074  Int_t ndim = GetNdimensions();
1075  TString name(GetName());
1076  for (Int_t d = 0; d < ndim; ++d)
1077  name += Form("_%d", group[d]);
1078 
1079 
1080  TString title(GetTitle());
1081  Ssiz_t posInsert = title.First(';');
1082  if (posInsert == kNPOS) {
1083  title += " rebin ";
1084  for (Int_t d = 0; d < ndim; ++d)
1085  title += Form("{%d}", group[d]);
1086  } else {
1087  for (Int_t d = ndim - 1; d >= 0; --d)
1088  title.Insert(posInsert, Form("{%d}", group[d]));
1089  title.Insert(posInsert, " rebin ");
1090  }
1091 
1092  TObjArray newaxes(ndim);
1093  newaxes.SetOwner();
1094  for (Int_t d = 0; d < ndim; ++d) {
1095  newaxes.AddAt(new TAxis(*GetAxis(d) ),d);
1096  if (group[d] > 1) {
1097  TAxis* newaxis = (TAxis*) newaxes.At(d);
1098  Int_t newbins = (newaxis->GetNbins() + group[d] - 1) / group[d];
1099  if (newaxis->GetXbins() && newaxis->GetXbins()->GetSize()) {
1100  // variable bins
1101  Double_t *edges = new Double_t[newbins + 1];
1102  for (Int_t i = 0; i < newbins + 1; ++i)
1103  if (group[d] * i <= newaxis->GetNbins())
1104  edges[i] = newaxis->GetXbins()->At(group[d] * i);
1105  else edges[i] = newaxis->GetXmax();
1106  newaxis->Set(newbins, edges);
1107  delete [] edges;
1108  } else {
1109  newaxis->Set(newbins, newaxis->GetXmin(), newaxis->GetXmax());
1110  }
1111  }
1112  }
1113 
1114  THnBase* h = CloneEmpty(name.Data(), title.Data(), &newaxes, kTRUE);
1115  Bool_t haveErrors = GetCalculateErrors();
1116  Bool_t wantErrors = haveErrors;
1117 
1118  Int_t* bins = new Int_t[ndim];
1119  Int_t* coord = new Int_t[fNdimensions];
1120 
1121  Long64_t i = 0;
1122  THnIter iter(this);
1123  while ((i = iter.Next(coord)) >= 0) {
1124  Double_t v = GetBinContent(i);
1125  for (Int_t d = 0; d < ndim; ++d) {
1126  bins[d] = TMath::CeilNint( (double) coord[d]/group[d] );
1127  }
1128  Long64_t idxh = h->GetBin(bins, kTRUE /*allocate*/);
1129 
1130  if (wantErrors) {
1131  // wantErrors == haveErrors, thus:
1132  h->AddBinError2(idxh, GetBinError2(i));
1133  }
1134 
1135  // only _after_ error calculation, or sqrt(v) is taken into account!
1136  h->AddBinContent(idxh, v);
1137  }
1138 
1139  delete [] bins;
1140  delete [] coord;
1141 
1142  h->SetEntries(fEntries);
1143 
1144  return h;
1145 
1146 }
1147 
1148 ////////////////////////////////////////////////////////////////////////////////
1149 /// Clear the histogram
1150 
1151 void THnBase::ResetBase(Option_t * /*option = ""*/)
1152 {
1153  fEntries = 0.;
1154  fTsumw = 0.;
1155  fTsumw2 = -1.;
1156  if (fIntegralStatus != kNoInt) {
1157  delete [] fIntegral;
1159  }
1160 }
1161 
1162 ////////////////////////////////////////////////////////////////////////////////
1163 /// Calculate the integral of the histogram
1164 
1166 {
1167  // delete old integral
1168  if (fIntegralStatus != kNoInt) {
1169  delete [] fIntegral;
1171  }
1172 
1173  // check number of bins
1174  if (GetNbins() == 0) {
1175  Error("ComputeIntegral", "The histogram must have at least one bin.");
1176  return 0.;
1177  }
1178 
1179  // allocate integral array
1180  fIntegral = new Double_t [GetNbins() + 1];
1181  fIntegral[0] = 0.;
1182 
1183  // fill integral array with contents of regular bins (non over/underflow)
1184  Int_t* coord = new Int_t[fNdimensions];
1185  Long64_t i = 0;
1186  THnIter iter(this);
1187  while ((i = iter.Next(coord)) >= 0) {
1188  Double_t v = GetBinContent(i);
1189 
1190  // check whether the bin is regular
1191  bool regularBin = true;
1192  for (Int_t dim = 0; dim < fNdimensions; dim++) {
1193  if (coord[dim] < 1 || coord[dim] > GetAxis(dim)->GetNbins()) {
1194  regularBin = false;
1195  break;
1196  }
1197  }
1198 
1199  // if outlayer, count it with zero weight
1200  if (!regularBin) v = 0.;
1201 
1202  fIntegral[i + 1] = fIntegral[i] + v;
1203  }
1204  delete [] coord;
1205 
1206  // check sum of weights
1207  if (fIntegral[GetNbins()] == 0.) {
1208  Error("ComputeIntegral", "No hits in regular bins (non over/underflow).");
1209  delete [] fIntegral;
1210  return 0.;
1211  }
1212 
1213  // normalize the integral array
1214  for (Long64_t j = 0; j <= GetNbins(); ++j)
1215  fIntegral[j] = fIntegral[j] / fIntegral[GetNbins()];
1216 
1217  // set status to valid
1219  return fIntegral[GetNbins()];
1220 }
1221 
1222 ////////////////////////////////////////////////////////////////////////////////
1223 /// Print bin with linex index "idx".
1224 /// For valid options see PrintBin(Long64_t idx, Int_t* bin, Option_t* options).
1225 
1226 void THnBase::PrintBin(Long64_t idx, Option_t* options) const
1227 {
1228  Int_t* coord = new Int_t[fNdimensions];
1229  PrintBin(idx, coord, options);
1230  delete [] coord;
1231 }
1232 
1233 ////////////////////////////////////////////////////////////////////////////////
1234 /// Print one bin. If "idx" is != -1 use that to determine the bin,
1235 /// otherwise (if "idx" == -1) use the coordinate in "bin".
1236 /// If "options" contains:
1237 /// '0': only print bins with an error or content != 0
1238 /// Return whether the bin was printed (depends on options)
1239 
1240 Bool_t THnBase::PrintBin(Long64_t idx, Int_t* bin, Option_t* options) const
1241 {
1242  Double_t v = -42;
1243  if (idx == -1) {
1244  idx = GetBin(bin);
1245  v = GetBinContent(idx);
1246  } else {
1247  v = GetBinContent(idx, bin);
1248  }
1249 
1250  Double_t err = 0.;
1251  if (GetCalculateErrors()) {
1252  if (idx != -1) {
1253  err = GetBinError(idx);
1254  }
1255  }
1256 
1257  if (v == 0. && err == 0. && options && strchr(options, '0')) {
1258  // suppress zeros, and we have one.
1259  return kFALSE;
1260  }
1261 
1262  TString coord;
1263  for (Int_t dim = 0; dim < fNdimensions; ++dim) {
1264  coord += bin[dim];
1265  coord += ',';
1266  }
1267  coord.Remove(coord.Length() - 1);
1268 
1269  if (GetCalculateErrors()) {
1270  Printf("Bin at (%s) = %g (+/- %g)", coord.Data(), v, err);
1271  } else {
1272  Printf("Bin at (%s) = %g", coord.Data(), v);
1273  }
1274 
1275  return kTRUE;
1276 }
1277 
1278 ////////////////////////////////////////////////////////////////////////////////
1279 /// Print "howmany" entries starting at "from". If "howmany" is -1, print all.
1280 /// If "options" contains:
1281 /// 'x': print in the order of axis bins, i.e. (0,0,...,0), (0,0,...,1),...
1282 /// '0': only print bins with content != 0
1283 
1284 void THnBase::PrintEntries(Long64_t from /*=0*/, Long64_t howmany /*=-1*/,
1285  Option_t* options /*=0*/) const
1286 {
1287  if (from < 0) from = 0;
1288  if (howmany == -1) howmany = GetNbins();
1289 
1290  Int_t* bin = new Int_t[fNdimensions];
1291 
1292  if (options && (strchr(options, 'x') || strchr(options, 'X'))) {
1293  Int_t* nbins = new Int_t[fNdimensions];
1294  for (Int_t dim = fNdimensions - 1; dim >= 0; --dim) {
1295  nbins[dim] = GetAxis(dim)->GetNbins();
1296  bin[dim] = from % nbins[dim];
1297  from /= nbins[dim];
1298  }
1299 
1300  for (Long64_t i = 0; i < howmany; ++i) {
1301  if (!PrintBin(-1, bin, options))
1302  ++howmany;
1303  // Advance to next bin:
1304  ++bin[fNdimensions - 1];
1305  for (Int_t dim = fNdimensions - 1; dim >= 0; --dim) {
1306  if (bin[dim] >= nbins[dim]) {
1307  bin[dim] = 0;
1308  if (dim > 0) {
1309  ++bin[dim - 1];
1310  } else {
1311  howmany = -1; // aka "global break"
1312  }
1313  }
1314  }
1315  }
1316  delete [] nbins;
1317  } else {
1318  for (Long64_t i = from; i < from + howmany; ++i) {
1319  if (!PrintBin(i, bin, options))
1320  ++howmany;
1321  }
1322  }
1323  delete [] bin;
1324 }
1325 
1326 ////////////////////////////////////////////////////////////////////////////////
1327 /// Print a THnBase. If "option" contains:
1328 /// 'a': print axis details
1329 /// 'm': print memory usage
1330 /// 's': print statistics
1331 /// 'c': print its content, too (this can generate a LOT of output!)
1332 /// Other options are forwarded to PrintEntries().
1333 
1334 void THnBase::Print(Option_t* options) const
1335 {
1336  Bool_t optAxis = options && (strchr(options, 'A') || (strchr(options, 'a')));
1337  Bool_t optMem = options && (strchr(options, 'M') || (strchr(options, 'm')));
1338  Bool_t optStat = options && (strchr(options, 'S') || (strchr(options, 's')));
1339  Bool_t optContent = options && (strchr(options, 'C') || (strchr(options, 'c')));
1340 
1341  Printf("%s (*0x%lx): \"%s\" \"%s\"", IsA()->GetName(), (unsigned long)this, GetName(), GetTitle());
1342  Printf(" %d dimensions, %g entries in %lld filled bins", GetNdimensions(), GetEntries(), GetNbins());
1343 
1344  if (optAxis) {
1345  for (Int_t dim = 0; dim < fNdimensions; ++dim) {
1346  TAxis* axis = GetAxis(dim);
1347  Printf(" axis %d \"%s\": %d bins (%g..%g), %s bin sizes", dim,
1348  axis->GetTitle(), axis->GetNbins(), axis->GetXmin(), axis->GetXmax(),
1349  (axis->GetXbins() ? "variable" : "fixed"));
1350  }
1351  }
1352 
1353  if (optStat) {
1354  Printf(" %s error calculation", (GetCalculateErrors() ? "with" : "without"));
1355  if (GetCalculateErrors()) {
1356  Printf(" Sum(w)=%g, Sum(w^2)=%g", GetSumw(), GetSumw2());
1357  for (Int_t dim = 0; dim < fNdimensions; ++dim) {
1358  Printf(" axis %d: Sum(w*x)=%g, Sum(w*x^2)=%g", dim, GetSumwx(dim), GetSumwx2(dim));
1359  }
1360  }
1361  }
1362 
1363  if (optMem && InheritsFrom(THnSparse::Class())) {
1364  const THnSparse* hsparse = dynamic_cast<const THnSparse*>(this);
1365  Printf(" coordinates stored in %d chunks of %d entries\n %g of bins filled using %g of memory compared to an array",
1366  hsparse->GetNChunks(), hsparse->GetChunkSize(),
1367  hsparse->GetSparseFractionBins(), hsparse->GetSparseFractionMem());
1368  }
1369 
1370  if (optContent) {
1371  Printf(" BIN CONTENT:");
1372  PrintEntries(0, -1, options);
1373  }
1374 }
1375 
1376 
1377 ////////////////////////////////////////////////////////////////////////////////
1378 /// Browse a THnSparse: create an entry (ROOT::THnSparseBrowsable) for each
1379 /// dimension.
1380 
1382 {
1383  if (fBrowsables.IsEmpty()) {
1384  for (Int_t dim = 0; dim < fNdimensions; ++dim) {
1386  }
1388  }
1389 
1390  for (Int_t dim = 0; dim < fNdimensions; ++dim) {
1391  b->Add(fBrowsables[dim]);
1392  }
1393 }
1394 
1395 
1396 
1397 
1398 /** \class ROOT::Internal::THnBaseBinIter
1399  Iterator over THnBase bins (internal implementation).
1400 */
1401 
1402 /// Destruct a bin iterator.
1403 
1405  // Not much to do, but pin vtable
1406 }
1407 
1408 
1409 
1410 /** \class THnIter
1411  Iterator over THnBase bins
1412 */
1413 
1414 ClassImp(THnIter);
1415 
1417  // Destruct a bin iterator.
1418  delete fIter;
1419 }
1420 
1421 
1422 /** \class ROOT::Internal::THnBaseBrowsable
1423  TBrowser helper for THnBase.
1424 */
1425 
1426 
1428 
1429 ////////////////////////////////////////////////////////////////////////////////
1430 /// Construct a THnBaseBrowsable.
1431 
1433 fHist(hist), fAxis(axis), fProj(0)
1434 {
1435  TString axisName = hist->GetAxis(axis)->GetName();
1436  if (axisName.IsNull()) {
1437  axisName = TString::Format("axis%d", axis);
1438  }
1439 
1440  SetNameTitle(axisName,
1441  TString::Format("Projection on %s of %s", axisName.Data(),
1442  hist->IsA()->GetName()).Data());
1443 }
1444 
1445 ////////////////////////////////////////////////////////////////////////////////
1446 /// Destruct a THnBaseBrowsable.
1447 
1449 {
1450  delete fProj;
1451 }
1452 
1453 ////////////////////////////////////////////////////////////////////////////////
1454 /// Browse an axis of a THnBase, i.e. draw its projection.
1455 
1457 {
1458  if (!fProj) {
1460  }
1461  fProj->Draw(b ? b->GetDrawOption() : "");
1462  gPad->Update();
1463 }
1464 
TString fTitle
Definition: TNamed.h:33
Abstract array base class.
Definition: TArray.h:31
void Add(TObject *obj, const char *name=0, Int_t check=-1)
Add object with name to browser.
Definition: TBrowser.cxx:261
Double_t ComputeIntegral()
Calculate the integral of the histogram.
Definition: THnBase.cxx:1165
Double_t GetBinError(const Int_t *idx) const
Definition: THnBase.h:159
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
Double_t GetSparseFractionMem() const
Return the amount of used memory over memory that would be used by a non-sparse n-dimensional histogr...
Definition: THnSparse.cxx:866
virtual Int_t GetNcells() const
Definition: TH1.h:280
virtual Double_t GetEffectiveEntries() const
Number of effective entries of the histogram.
Definition: TH1.cxx:4079
An array of TObjects.
Definition: TObjArray.h:37
float xmin
Definition: THbookFile.cxx:93
virtual void Reset(Option_t *option="")=0
Bool_t GetCalculateErrors() const
Definition: THnBase.h:136
void Divide(const THnBase *h)
Divide this histogram by h this = this/(h) Note that if h has Sumw2 set, Sumw2 is automatically calle...
Definition: THnBase.cxx:847
Double_t Floor(Double_t x)
Definition: TMath.h:600
Int_t GetFirst() const
Return first bin on the axis i.e.
Definition: TAxis.cxx:444
long long Long64_t
Definition: RtypesCore.h:69
Double_t * fIntegral
Definition: THnBase.h:53
void SetEntries(Double_t entries)
Definition: THnBase.h:164
TFitResultPtr Fit(TF1 *f1, Option_t *option="", Option_t *goption="")
Fit a THnSparse with function f.
Definition: THnBase.cxx:365
Double_t fTsumw2
Definition: THnBase.h:50
THnBase * RebinBase(Int_t group) const
Combine the content of "group" neighboring bins into a new bin and return the resulting THnBase...
Definition: THnBase.cxx:1056
virtual Long64_t GetNbins() const =0
const char Option_t
Definition: RtypesCore.h:62
void Multiply(const THnBase *h)
Multiply this histogram by histogram h this = this * h Note that if h has Sumw2 set, Sumw2 is automatically called for this if not already set.
Definition: THnBase.cxx:756
virtual Double_t GetBinLowEdge(Int_t bin) const
Return low edge of bin.
Definition: TAxis.cxx:504
return c1
Definition: legend1.C:41
const Ssiz_t kNPOS
Definition: RtypesCore.h:115
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:640
const Double_t * GetArray() const
Definition: TArrayD.h:43
Bool_t TestBit(UInt_t f) const
Definition: TObject.h:159
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition: TNamed.cxx:131
virtual void SetOwner(Bool_t enable=kTRUE)
Set whether this collection is the owner (enable==true) of its content.
TH1 * h
Definition: legend2.C:5
virtual ~THnIter()
Definition: THnBase.cxx:1416
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH1.cxx:4639
Double_t GetBinContent(const Int_t *idx) const
Definition: THnBase.h:167
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition: TString.h:587
void Scale(Double_t c)
Scale contents and errors of this histogram by c: this = this * c It does not modify the histogram&#39;s ...
Definition: THnBase.cxx:610
virtual void GetBinXYZ(Int_t binglobal, Int_t &binx, Int_t &biny, Int_t &binz) const
Return binx, biny, binz corresponding to the global bin number globalbin see TH1::GetBin function abo...
Definition: TH1.cxx:4554
static THnBase * CreateHnAny(const char *name, const char *title, const TH1 *h1, Bool_t sparse, Int_t chunkSize=1024 *16)
Create a THn / THnSparse object from a histogram deriving from TH1.
Definition: THnBase.cxx:197
void GetRandom(Double_t *rand, Bool_t subBinRandom=kTRUE)
Generate an n-dimensional random tuple based on the histogrammed distribution.
Definition: THnBase.cxx:394
Basic string class.
Definition: TString.h:129
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1099
Bool_t IsEmpty() const
Definition: TObjArray.h:70
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
int nbins[3]
#define NULL
Definition: RtypesCore.h:88
TObject * At(Int_t idx) const
Definition: TObjArray.h:165
Int_t GetNChunks() const
Definition: THnSparse.h:87
virtual void SetNameTitle(const char *name, const char *title)
Set all the TNamed parameters (name and title).
Definition: TNamed.cxx:145
TString & Insert(Ssiz_t pos, const char *s)
Definition: TString.h:597
void AddBinContent(const Int_t *x, Double_t v=1.)
Definition: THnBase.h:163
Short_t Abs(Short_t d)
Definition: TMathBase.h:108
TObjArray fAxes
Definition: THnBase.h:46
int Nostore
Definition: Foption.h:41
void Reset()
Definition: TCollection.h:156
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:687
Double_t fTsumw
Definition: THnBase.h:49
virtual Double_t GetBinUpEdge(Int_t bin) const
Return up edge of bin.
Definition: TAxis.cxx:514
virtual Int_t GetDimension() const
Definition: TH1.h:263
Double_t fEntries
browser-helpers for each axis
Definition: THnBase.h:48
THnBaseBrowsable(THnBase *hist, Int_t axis)
Construct a THnBaseBrowsable.
Definition: THnBase.cxx:1432
Double_t GetXmin() const
Definition: TAxis.h:133
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition: TObject.cxx:135
Double_t x[n]
Definition: legend1.C:17
TBrowser helper for THnBase.
Definition: THnBase.h:271
static Int_t FitOptionsMake(Option_t *option, Foption_t &Foption)
Decode string choptin and fill fitOption structure.
Definition: TH1.cxx:4232
Efficient multidimensional histogram.
Definition: THnSparse.h:36
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString...
Definition: TString.cxx:2345
void Class()
Definition: Class.C:29
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:29
TH1D * Projection(Int_t xDim, Option_t *option="") const
Definition: THnBase.h:188
int Chi2
Definition: Foption.h:32
THnBase * CloneEmpty(const char *name, const char *title, const TObjArray *axes, Bool_t keepTargetAxis) const
Create a new THnBase object that is of the same type as *this, but with dimensions and bins given by ...
Definition: THnBase.cxx:80
Int_t GetCoord(Int_t dim) const
Definition: THnBase.h:318
virtual void Reserve(Long64_t)
Definition: THnBase.h:99
virtual Double_t GetBinCenter(Int_t bin) const
Return center of bin.
Definition: TAxis.cxx:464
void PrintEntries(Long64_t from=0, Long64_t howmany=-1, Option_t *options=0) const
Print "howmany" entries starting at "from".
Definition: THnBase.cxx:1284
void SetBinContent(const Int_t *idx, Double_t v)
Definition: THnBase.h:177
std::vector< std::vector< double > > Data
Double_t GetSumw2() const
Definition: THnBase.h:184
virtual ~THnBaseBinIter()
Destruct a bin iterator.
Definition: THnBase.cxx:1404
Double_t GetSumw() const
Definition: THnBase.h:183
TObjArray fBrowsables
Definition: THnBase.h:47
Ssiz_t First(char c) const
Find first occurrence of a character c.
Definition: TString.cxx:477
Bool_t CheckConsistency(const THnBase *h, const char *tag) const
Consistency check on (some of) the parameters of two histograms (for operations). ...
Definition: THnBase.cxx:986
void Browse(TBrowser *b)
Browse an axis of a THnBase, i.e. draw its projection.
Definition: THnBase.cxx:1456
TH1F * h1
Definition: legend1.C:5
virtual void AddBinContent(Int_t bin)
Increment bin content by 1.
Definition: TH1.cxx:1193
virtual void ResetStats()
Reset the statistics including the number of entries and replace with values calculates from bin cont...
Definition: TH1.cxx:7088
virtual void SetBinError(Int_t bin, Double_t error)
See convention for numbering bins in TH1::GetBin.
Definition: TH1.cxx:8311
virtual Double_t GetBinError2(Long64_t linidx) const =0
Bool_t PrintBin(Long64_t idx, Int_t *coord, Option_t *options) const
Print one bin.
Definition: THnBase.cxx:1240
virtual Double_t Rndm()
Machine independent random number generator.
Definition: TRandom.cxx:512
THist< 3, double, THistStatContent, THistStatUncertainty > TH3D
Definition: THist.hxx:322
Using a TBrowser one can browse all ROOT objects.
Definition: TBrowser.h:37
RooArgSet S(const RooAbsArg &v1)
void AddInternal(const THnBase *h, Double_t c, Bool_t rebinned)
Add() implementation for both rebinned histograms and those with identical binning.
Definition: THnBase.cxx:634
#define F(x, y, z)
virtual void SetRange(Int_t first=0, Int_t last=0)
Set the viewing range for the axis from bin first to last.
Definition: TAxis.cxx:889
virtual void InitStorage(Int_t *nbins, Int_t chunkSize)=0
virtual void SetFilledBins(Long64_t)
Definition: THnBase.h:100
static double C[]
TAxis * GetAxis(Int_t dim) const
Definition: THnBase.h:125
#define R__THNBCASE(TAG)
virtual void AddAtAndExpand(TObject *obj, Int_t idx)
Add object at position idx.
Definition: TObjArray.cxx:222
void Print(Option_t *option="") const
Print a THnBase.
Definition: THnBase.cxx:1334
const char * GetTitle() const
Returns title of object.
Definition: TAxis.h:129
Int_t GetLast() const
Return last bin on the axis i.e.
Definition: TAxis.cxx:455
Class to manage histogram axis.
Definition: TAxis.h:30
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2851
Int_t GetSize() const
Definition: TArray.h:47
SVector< double, 2 > v
Definition: Dict.h:5
Long64_t Next(Int_t *coord=0)
Definition: THnBase.h:310
virtual Long64_t GetBin(const Int_t *idx) const =0
void Init(const char *name, const char *title, const TObjArray *axes, Bool_t keepTargetAxis, Int_t chunkSize=1024 *16)
Initialize axes and name.
Definition: THnBase.cxx:96
Provides an indirection to the TFitResult class and with a semantics identical to a TFitResult pointe...
Definition: TFitResultPtr.h:31
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:436
Collection abstract base class.
Definition: TCollection.h:42
virtual ~THnBase()
Destruct a THnBase.
Definition: THnBase.cxx:69
Int_t GetEntriesFast() const
Definition: TObjArray.h:64
Double_t GetSumwx2(Int_t dim) const
Definition: THnBase.h:186
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:873
char * Form(const char *fmt,...)
Ssiz_t Length() const
Definition: TString.h:388
static void RejectPoint(Bool_t reject=kTRUE)
Static function to set the global flag to reject points the fgRejectPoint global flag is tested by al...
Definition: TF1.cxx:3459
int Like
Definition: Foption.h:34
The ROOT global object gROOT contains a list of all defined classes.
Definition: TClass.h:71
TAxis * GetYaxis()
Definition: TH1.h:301
double maxRange[3]
float xmax
Definition: THbookFile.cxx:93
void SetBinEdges(Int_t idim, const Double_t *bins)
Set the axis # of bins and bin limits on dimension idim.
Definition: THnBase.cxx:1004
R__EXTERN TRandom * gRandom
Definition: TRandom.h:62
virtual Int_t GetBin(Int_t binx, Int_t biny=0, Int_t binz=0) const
Return Global bin number corresponding to binx,y,z.
Definition: TH1.cxx:4541
virtual void Rebuild(Option_t *option="")
Using the current bin info, recompute the arrays for contents and errors.
Definition: TH1.cxx:6403
virtual void AddAt(TObject *obj, Int_t idx)
Add object at position ids.
Definition: TObjArray.cxx:239
Double_t GetEntries() const
Definition: THnBase.h:133
#define Printf
Definition: TGeoToOCC.h:18
const Bool_t kFALSE
Definition: RtypesCore.h:92
Double_t At(Int_t i) const
Definition: TArrayD.h:79
Double_t GetSumwx(Int_t dim) const
Definition: THnBase.h:185
TString & Remove(Ssiz_t pos)
Definition: TString.h:621
void AddRange(unsigned int icoord, double xmin, double xmax)
add a range [xmin,xmax] for the new coordinate icoord Adding a range does not delete existing one...
Definition: DataRange.cxx:94
int Ssiz_t
Definition: RtypesCore.h:63
void SetBinError(const Int_t *idx, Double_t e)
Definition: THnBase.h:161
virtual Bool_t IsEmpty() const
Definition: TCollection.h:93
return c2
Definition: legend2.C:14
virtual void SetBinError2(Long64_t bin, Double_t e2)=0
#define ClassImp(name)
Definition: Rtypes.h:336
class describing the range in the coordinates it supports multiple range in a coordinate.
Definition: DataRange.h:34
double f(double x)
double Double_t
Definition: RtypesCore.h:55
Bool_t IsInRange(Int_t *coord) const
Check whether bin coord is in range, as defined by TAxis::SetRange().
Definition: THnBase.cxx:429
void Add(const THnBase *h, Double_t c=1.)
Add contents of h scaled by c to this histogram: this = this + c * h Note that if h has Sumw2 set...
Definition: THnBase.cxx:695
Long64_t Merge(TCollection *list)
Merge this with a list of THnBase&#39;s.
Definition: THnBase.cxx:721
int type
Definition: TGX11.cxx:120
THnBase()
Definition: THnBase.h:65
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:572
The TH1 histogram class.
Definition: TH1.h:56
~THnBaseBrowsable()
Destruct a THnBaseBrowsable.
Definition: THnBase.cxx:1448
enum THnBase::@97 fIntegralStatus
array with bin weight sums
THist< 2, double, THistStatContent, THistStatUncertainty > TH2D
Definition: THist.hxx:316
virtual void AddBinError2(Long64_t bin, Double_t e2)=0
Double_t GetSparseFractionBins() const
Return the amount of filled bins over all bins.
Definition: THnSparse.cxx:855
virtual void Sumw2()=0
static TClass * GetClass(const char *name, Bool_t load=kTRUE, Bool_t silent=kFALSE)
Static method returning pointer to TClass of the specified class name.
Definition: TClass.cxx:2885
TObject * ProjectionAny(Int_t ndim, const Int_t *dim, Bool_t wantNDim, Option_t *option="") const
Project all bins into a ndim-dimensional THn / THnSparse (whatever *this is) or if (ndim < 4 and !wan...
Definition: THnBase.cxx:454
Bool_t IsNull() const
Definition: TString.h:385
Int_t GetNdimensions() const
Definition: THnBase.h:135
TAxis * GetZaxis()
Definition: TH1.h:302
void Browse(TBrowser *b)
Browse a THnSparse: create an entry (ROOT::THnSparseBrowsable) for each dimension.
Definition: THnBase.cxx:1381
Mother of all ROOT objects.
Definition: TObject.h:37
virtual Bool_t IsInside(const Double_t *x) const
return kTRUE if the point is inside the function range
Definition: TF1.h:542
virtual Double_t GetBinWidth(Int_t bin) const
Return bin width.
Definition: TAxis.cxx:526
Int_t fNdimensions
Definition: THnBase.h:45
Bool_t HaveSkippedBin() const
Definition: THnBase.h:319
Iterator over THnBase bins.
Definition: THnBase.h:304
1-Dim function class
Definition: TF1.h:150
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
THist< 1, double, THistStatContent, THistStatUncertainty > TH1D
Definition: THist.hxx:310
#define gPad
Definition: TVirtualPad.h:284
TObjArray * GetListOfAxes()
Definition: THnBase.h:123
virtual void SetEntries(Double_t n)
Definition: TH1.h:363
Multidimensional histogram base.
Definition: THnBase.h:43
void ResetBit(UInt_t f)
Definition: TObject.h:158
TFitResultPtr FitObject(TH1 *h1, TF1 *f1, Foption_t &option, const ROOT::Math::MinimizerOptions &moption, const char *goption, ROOT::Fit::DataRange &range)
fitting function for a TH1 (called from TH1::Fit)
Definition: HFitImpl.cxx:943
Option_t * GetDrawOption() const
Get option used by the graphics system to draw this object.
Definition: TBrowser.h:104
Double_t Sqrt(Double_t x)
Definition: TMath.h:591
Int_t GetNbins() const
Definition: TAxis.h:121
#define I(x, y, z)
TH1 * CreateHist(const char *name, const char *title, const TObjArray *axes, Bool_t keepTargetAxis) const
Create an empty histogram with name and title with a given set of axes.
Definition: THnBase.cxx:144
virtual const char * GetName() const
Returns name of object.
Definition: TObject.cxx:364
virtual Int_t GetSize() const
Definition: TCollection.h:89
virtual void Set(Int_t nbins, Double_t xmin, Double_t xmax)
Initialize axis with fix bins.
Definition: TAxis.cxx:717
Int_t GetChunkSize() const
Definition: THnSparse.h:86
virtual Double_t EvalPar(const Double_t *x, const Double_t *params=0)
Evaluate function with given coordinates and parameters.
Definition: TF1.cxx:1226
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition: TNamed.cxx:155
const Bool_t kTRUE
Definition: RtypesCore.h:91
Double_t GetXmax() const
Definition: TAxis.h:134
void RebinnedAdd(const THnBase *h, Double_t c=1.)
Add contents of h scaled by c to this histogram: this = this + c * h Note that if h has Sumw2 set...
Definition: THnBase.cxx:711
double minRange[3]
Int_t CeilNint(Double_t x)
Definition: TMath.h:597
Long64_t BinarySearch(Long64_t n, const T *array, T value)
Definition: TMath.h:1093
const TArrayD * GetXbins() const
Definition: TAxis.h:130
void ResetBase(Option_t *option="")
Clear the histogram.
Definition: THnBase.cxx:1151
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:859
TAxis * GetXaxis()
Definition: TH1.h:300
void SetTitle(const char *title)
Change (i.e.
Definition: THnBase.cxx:1020
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
virtual Double_t GetBinError(Int_t bin) const
Return value of error associated to bin number bin.
Definition: TH1.cxx:8175
void * New(ENewType defConstructor=kClassNew, Bool_t quiet=kFALSE) const
Return a pointer to a newly allocated object of this class.
Definition: TClass.cxx:4706
const char * Data() const
Definition: TString.h:347