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