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