Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
THnBase.h
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#ifndef ROOT_THnBase
13#define ROOT_THnBase
14
15/*************************************************************************
16
17 THnBase: Common base class for n-dimensional histogramming.
18 Defines interfaces and algorithms.
19
20*************************************************************************/
21
22
23#include "TNamed.h"
24#include "TMath.h"
25#include "TFitResultPtr.h"
26#include "TObjArray.h"
27#include "TArrayD.h"
28
29class TAxis;
30class TH1;
31class TH1D;
32class TH2D;
33class TH3D;
34class TF1;
35class THnIter;
36
37namespace ROOT {
38namespace Internal {
39 class THnBaseBinIter;
40}
41}
42
43class THnBase: public TNamed {
44protected:
45 Int_t fNdimensions; ///< Number of dimensions
46 TObjArray fAxes; ///< Axes of the histogram
47 TObjArray fBrowsables; ///<! Browser-helpers for each axis
48 Double_t fEntries; ///< Number of entries, spread over chunks
49 Double_t fTsumw; ///< Total sum of weights
50 Double_t fTsumw2; ///< Total sum of weights squared; -1 if no errors are calculated
51 TArrayD fTsumwx; ///< Total sum of weight*X for each dimension
52 TArrayD fTsumwx2; ///< Total sum of weight*X*X for each dimension
53 std::vector<Double_t> fIntegral; ///<! vector with bin weight sums
54 enum {
58 } fIntegralStatus; ///<! status of integral
59
60 protected:
62
63 THnBase(const char *name, const char *title, Int_t dim, const Int_t *nbins, const Double_t *xmin,
64 const Double_t *xmax);
65
66 THnBase(const char *name, const char *title, Int_t dim, const Int_t *nbins,
67 const std::vector<std::vector<double>> &xbins);
68
69 THnBase(const THnBase &other);
70
71 THnBase &operator=(const THnBase &other);
72
73 THnBase(THnBase &&other);
74
75 THnBase &operator=(THnBase &&other);
76
77 void UpdateXStat(const Double_t *x, Double_t w = 1.)
78 {
79 if (GetCalculateErrors()) {
80 for (Int_t d = 0; d < fNdimensions; ++d) {
81 const Double_t xd = x[d];
82 fTsumwx[d] += w * xd;
83 fTsumwx2[d] += w * xd * xd;
84 }
85 }
86 }
87
88 /// Increment the statistics due to filled weight "w",
90 fEntries += 1;
91 if (GetCalculateErrors()) {
92 fTsumw += w;
93 fTsumw2 += w*w;
94 }
96 }
97
98 virtual void InitStorage(Int_t* nbins, Int_t chunkSize) = 0;
99 void Init(const char* name, const char* title,
100 const TObjArray* axes, Bool_t keepTargetAxis,
101 Int_t chunkSize = 1024 * 16);
102 THnBase* CloneEmpty(const char* name, const char* title,
103 const TObjArray* axes, Bool_t keepTargetAxis) const;
104 virtual void Reserve(Long64_t /*nbins*/) {}
105 virtual void SetFilledBins(Long64_t /*nbins*/) {};
106
107 Bool_t CheckConsistency(const THnBase *h, const char *tag) const;
108 TH1* CreateHist(const char* name, const char* title,
109 const TObjArray* axes, Bool_t keepTargetAxis) const;
110 TObject* ProjectionAny(Int_t ndim, const Int_t* dim,
111 Bool_t wantNDim, Option_t* option = "") const;
112 Bool_t PrintBin(Long64_t idx, Int_t* coord, Option_t* options) const;
113 void AddInternal(const THnBase* h, Double_t c, Bool_t rebinned);
115 THnBase* RebinBase(const Int_t* group) const;
116 void ResetBase(Option_t *option= "");
117
118 static THnBase* CreateHnAny(const char* name, const char* title,
119 const TH1* h1, Bool_t sparse,
120 Int_t chunkSize = 1024 * 16);
121 static THnBase* CreateHnAny(const char* name, const char* title,
122 const THnBase* hn, Bool_t sparse,
123 Int_t chunkSize = 1024 * 16);
124
125 public:
126 ~THnBase() override;
127
129 const TObjArray* GetListOfAxes() const { return &fAxes; }
130 TAxis* GetAxis(Int_t dim) const { return (TAxis*)fAxes[dim]; }
131
132 TFitResultPtr Fit(TF1 *f1 ,Option_t *option = "", Option_t *goption = "");
133 TList* GetListOfFunctions() { return nullptr; }
134
135 virtual ROOT::Internal::THnBaseBinIter* CreateIter(Bool_t respectAxisRange) const = 0;
136
137 virtual Long64_t GetNbins() const = 0;
138 Double_t GetEntries() const { return fEntries; }
139 Double_t GetWeightSum() const { return fTsumw; }
141 Bool_t GetCalculateErrors() const { return fTsumw2 >= 0.; }
142
143 /// Calculate errors (or not if "calc" == kFALSE)
145 if (calc) Sumw2();
146 else fTsumw2 = -1.;
147 }
148
150 UpdateXStat(x, w);
151 Long64_t bin = GetBin(x, kTRUE /*alloc*/);
152 FillBin(bin, w);
153 return bin;
154 }
155 Long64_t Fill(const char* name[], Double_t w = 1.) {
156 Long64_t bin = GetBin(name, kTRUE /*alloc*/);
157 FillBin(bin, w);
158 return bin;
159 }
160
161 /// Fill with the provided variadic arguments.
162 /// The number of arguments must be equal to the number of histogram dimensions or, for weighted fills, to the
163 /// number of dimensions + 1; in the latter case, the last function argument is used as weight.
164 /// A separate `firstval` argument is needed so the compiler does not pick this overload instead of the non-templated
165 /// Fill overloads
166 template <typename... MoreTypes>
167 Long64_t Fill(Double_t firstval, MoreTypes... morevals)
168 {
169 const std::array<double, 1 + sizeof...(morevals)> x{firstval, static_cast<double>(morevals)...};
170 if (Int_t(x.size()) == GetNdimensions()) {
171 // without weight
172 return Fill(x.data());
173 } else if (Int_t(x.size()) == (GetNdimensions() + 1)) {
174 // with weight
175 return Fill(x.data(), x.back());
176 } else {
177 Error("Fill", "Wrong number of arguments for number of histogram axes.");
178 }
179
180 return -1;
181 }
182
183 virtual void FillBin(Long64_t bin, Double_t w) = 0;
184
185 void SetBinEdges(Int_t idim, const Double_t* bins);
186 Bool_t IsInRange(Int_t *coord) const;
187 Double_t GetBinError(const Int_t *idx) const { return GetBinError(GetBin(idx)); }
188 Double_t GetBinError(Long64_t linidx) const { return TMath::Sqrt(GetBinError2(linidx)); }
189 void SetBinError(const Int_t* idx, Double_t e) { SetBinError(GetBin(idx), e); }
191 void AddBinContent(const Int_t* x, Double_t v = 1.) { AddBinContent(GetBin(x), v); }
192 void SetEntries(Double_t entries) { fEntries = entries; }
193 void SetTitle(const char *title) override;
194
195 Double_t GetBinContent(const Int_t *idx) const { return GetBinContent(GetBin(idx)); } // intentionally non-virtual
196 virtual Double_t GetBinContent(Long64_t bin, Int_t* idx = nullptr) const = 0;
197 virtual Double_t GetBinError2(Long64_t linidx) const = 0;
198 virtual Long64_t GetBin(const Int_t* idx) const = 0;
199 virtual Long64_t GetBin(const Double_t* x) const = 0;
200 virtual Long64_t GetBin(const char* name[]) const = 0;
201 virtual Long64_t GetBin(const Int_t* idx, Bool_t /*allocate*/ = kTRUE) = 0;
202 virtual Long64_t GetBin(const Double_t* x, Bool_t /*allocate*/ = kTRUE) = 0;
203 virtual Long64_t GetBin(const char* name[], Bool_t /*allocate*/ = kTRUE) = 0;
204
205 void SetBinContent(const Int_t* idx, Double_t v) { SetBinContent(GetBin(idx), v); } // intentionally non-virtual
206 virtual void SetBinContent(Long64_t bin, Double_t v) = 0;
207 virtual void SetBinError2(Long64_t bin, Double_t e2) = 0;
208 virtual void AddBinError2(Long64_t bin, Double_t e2) = 0;
209 virtual void AddBinContent(Long64_t bin, Double_t v = 1.) = 0;
210
211 Double_t GetSumw() const { return fTsumw; }
212 Double_t GetSumw2() const { return fTsumw2; }
213 Double_t GetSumwx(Int_t dim) const { return fTsumwx[dim]; }
214 Double_t GetSumwx2(Int_t dim) const { return fTsumwx2[dim]; }
215
216 /// Project all bins into a 1-dimensional histogram,
217 /// keeping only axis "xDim".
218 /// If "option" contains:
219 /// - "E" errors will be calculated.
220 /// - "A" ranges of the taget axes will be ignored.
221 /// - "O" original axis range of the taget axes will be
222 /// kept, but only bins inside the selected range
223 /// will be filled.
224 TH1D* Projection(Int_t xDim, Option_t* option = "") const {
225 return (TH1D*) ProjectionAny(1, &xDim, false, option);
226 }
227
228 /// Project all bins into a 2-dimensional histogram,
229 /// keeping only axes "xDim" and "yDim".
230 ///
231 /// WARNING: just like TH3::Project3D("yx") and TTree::Draw("y:x"),
232 /// Projection(y,x) uses the first argument to define the y-axis and the
233 /// second for the x-axis!
234 ///
235 /// If "option" contains "E" errors will be calculated.
236 /// "A" ranges of the taget axes will be ignored.
237 TH2D* Projection(Int_t yDim, Int_t xDim, Option_t* option = "") const {
238 const Int_t dim[2] = {xDim, yDim};
239 return (TH2D*) ProjectionAny(2, dim, false, option);
240 }
241
242 /// Project all bins into a 3-dimensional histogram,
243 /// keeping only axes "xDim", "yDim", and "zDim".
244 /// If "option" contains:
245 /// - "E" errors will be calculated.
246 /// - "A" ranges of the taget axes will be ignored.
247 /// - "O" original axis range of the taget axes will be
248 /// kept, but only bins inside the selected range
249 /// will be filled.
250 TH3D* Projection(Int_t xDim, Int_t yDim, Int_t zDim, Option_t* option = "") const {
251 const Int_t dim[3] = {xDim, yDim, zDim};
252 return (TH3D*) ProjectionAny(3, dim, false, option);
253 }
254
255 THnBase* ProjectionND(Int_t ndim, const Int_t* dim,
256 Option_t* option = "") const {
257 return (THnBase*)ProjectionAny(ndim, dim, kTRUE /*wantNDim*/, option);
258 }
259
261
262 void Scale(Double_t c);
263 void Add(const THnBase* h, Double_t c=1.);
264 void Add(const TH1* hist, Double_t c=1.);
265 void Multiply(const THnBase* h);
266 void Multiply(TF1* f, Double_t c = 1.);
267 void Divide(const THnBase* h);
268 void Divide(const THnBase* h1, const THnBase* h2, Double_t c1 = 1., Double_t c2 = 1., Option_t* option="");
269 void RebinnedAdd(const THnBase* h, Double_t c=1.);
270
271 virtual void Reset(Option_t* option = "") = 0;
272 virtual void Sumw2() = 0;
273
275 void GetRandom(Double_t *rand, Bool_t subBinRandom = kTRUE);
276
277 void Print(Option_t* option = "") const override;
278 void PrintEntries(Long64_t from = 0, Long64_t howmany = -1, Option_t* options = nullptr) const;
279 void PrintBin(Int_t* coord, Option_t* options) const {
280 PrintBin(-1, coord, options);
281 }
282 void PrintBin(Long64_t idx, Option_t* options) const;
283
284 void Browse(TBrowser *b) override;
285 Bool_t IsFolder() const override { return kTRUE; }
286
287 //void Draw(Option_t* option = "");
288
289 ClassDefOverride(THnBase, 1); // Common base for n-dimensional histogram
290
291 friend class THnIter;
292};
293
294namespace ROOT {
295namespace Internal {
296 // Helper class for browsing THnBase objects
297 class THnBaseBrowsable: public TNamed {
298 public:
299 THnBaseBrowsable(THnBase* hist, Int_t axis);
300 ~THnBaseBrowsable() override;
301 void Browse(TBrowser *b) override;
302 Bool_t IsFolder() const override { return kFALSE; }
303
304 private:
305 THnBase* fHist; // Original histogram
306 Int_t fAxis; // Axis to visualize
307 TH1* fProj; // Projection result
308 ClassDefOverride(THnBaseBrowsable, 0); // Browser-helper for THnBase
309 };
310
311 // Base class for iterating over THnBase bins
313 public:
314 THnBaseBinIter(Bool_t respectAxisRange):
315 fRespectAxisRange(respectAxisRange), fHaveSkippedBin(kFALSE) {}
316 virtual ~THnBaseBinIter();
319
320 virtual Int_t GetCoord(Int_t dim) const = 0;
321 virtual Long64_t Next(Int_t* coord = nullptr) = 0;
322
323 protected:
326 };
327}
328}
329
330class THnIter: public TObject {
331public:
332 THnIter(const THnBase* hist, Bool_t respectAxisRange = kFALSE):
333 fIter(hist->CreateIter(respectAxisRange)) {}
334 ~THnIter() override;
335
336 /// Return the next bin's index.
337 /// If provided, set coord to that bin's coordinates (bin indexes).
338 /// I.e. coord must point to Int_t[hist->GetNdimensions()]
339 /// Returns -1 when all bins have been visited.
340 Long64_t Next(Int_t* coord = nullptr) {
341 return fIter->Next(coord);
342 }
343
344 Int_t GetCoord(Int_t dim) const { return fIter->GetCoord(dim); }
347
348private:
350 ClassDefOverride(THnIter, 0); //Iterator over bins of a THnBase.
351};
352
353#endif // ROOT_THnBase
#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
#define e(i)
Definition RSha256.hxx:103
bool Bool_t
Definition RtypesCore.h:63
int Int_t
Definition RtypesCore.h:45
constexpr Bool_t kFALSE
Definition RtypesCore.h:101
double Double_t
Definition RtypesCore.h:59
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 ClassDefOverride(name, id)
Definition Rtypes.h:341
Option_t Option_t option
char name[80]
Definition TGX11.cxx:110
float xmin
float xmax
Iterator over THnBase bins (internal implementation).
Definition THnBase.h:312
Bool_t RespectsAxisRange() const
Definition THnBase.h:318
Bool_t HaveSkippedBin() const
Definition THnBase.h:317
virtual Int_t GetCoord(Int_t dim) const =0
THnBaseBinIter(Bool_t respectAxisRange)
Definition THnBase.h:314
virtual ~THnBaseBinIter()
Destruct a bin iterator.
Definition THnBase.cxx:1516
virtual Long64_t Next(Int_t *coord=nullptr)=0
TBrowser helper for THnBase.
Definition THnBase.h:297
~THnBaseBrowsable() override
Destruct a THnBaseBrowsable.
Definition THnBase.cxx:1560
Bool_t IsFolder() const override
Returns kTRUE in case object contains browsable objects (like containers or lists of other objects).
Definition THnBase.h:302
void Browse(TBrowser *b) override
Browse an axis of a THnBase, i.e. draw its projection.
Definition THnBase.cxx:1568
Array of doubles (64 bits per element).
Definition TArrayD.h:27
Class to manage histogram axis.
Definition TAxis.h:30
Using a TBrowser one can browse all ROOT objects.
Definition TBrowser.h:37
Collection abstract base class.
Definition TCollection.h:65
1-Dim function class
Definition TF1.h:213
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:620
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:58
2-D histogram with a double per channel (see TH1 documentation)}
Definition TH2.h:300
3-D histogram with a double per channel (see TH1 documentation)}
Definition TH3.h:307
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:1493
const TObjArray * GetListOfAxes() const
Definition THnBase.h:129
Long64_t Fill(const Double_t *x, Double_t w=1.)
Definition THnBase.h:149
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
enum THnBase::@84 fIntegralStatus
! status of integral
Double_t fEntries
Number of entries, spread over chunks.
Definition THnBase.h:48
virtual void AddBinContent(Long64_t bin, Double_t v=1.)=0
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:528
TList * GetListOfFunctions()
Definition THnBase.h:133
virtual Long64_t GetBin(const Int_t *idx, Bool_t=kTRUE)=0
TFitResultPtr Fit(TF1 *f1, Option_t *option="", Option_t *goption="")
Fit a THnSparse with function f.
Definition THnBase.cxx:464
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:710
virtual void SetBinError2(Long64_t bin, Double_t e2)=0
TObjArray * GetListOfAxes()
Definition THnBase.h:128
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:243
Bool_t PrintBin(Long64_t idx, Int_t *coord, Option_t *options) const
Print one bin.
Definition THnBase.cxx:1352
TH1D * Projection(Int_t xDim, Option_t *option="") const
Project all bins into a 1-dimensional histogram, keeping only axis "xDim".
Definition THnBase.h:224
Double_t GetSumwx(Int_t dim) const
Definition THnBase.h:213
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:1098
Int_t GetNdimensions() const
Definition THnBase.h:140
TH3D * Projection(Int_t xDim, Int_t yDim, Int_t zDim, Option_t *option="") const
Project all bins into a 3-dimensional histogram, keeping only axes "xDim", "yDim",...
Definition THnBase.h:250
~THnBase() override
Destruct a THnBase.
Definition THnBase.cxx:165
void PrintBin(Int_t *coord, Option_t *options) const
Definition THnBase.h:279
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:296
virtual Long64_t GetBin(const char *name[], Bool_t=kTRUE)=0
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:1263
virtual Long64_t GetNbins() const =0
Bool_t IsFolder() const override
Returns kTRUE in case object contains browsable objects (like containers or lists of other objects).
Definition THnBase.h:285
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:1168
THnBase * ProjectionND(Int_t ndim, const Int_t *dim, Option_t *option="") const
Definition THnBase.h:255
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:1116
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:734
void SetTitle(const char *title) override
Change (i.e.
Definition THnBase.cxx:1132
void PrintEntries(Long64_t from=0, Long64_t howmany=-1, Option_t *options=nullptr) const
Print "howmany" entries starting at "from".
Definition THnBase.cxx:1396
virtual void SetBinContent(Long64_t bin, Double_t v)=0
Double_t GetBinError(Long64_t linidx) const
Definition THnBase.h:188
Double_t GetBinContent(const Int_t *idx) const
Definition THnBase.h:195
virtual Long64_t GetBin(const Double_t *x, Bool_t=kTRUE)=0
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:807
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
Long64_t Fill(Double_t firstval, MoreTypes... morevals)
Fill with the provided variadic arguments.
Definition THnBase.h:167
virtual Long64_t GetBin(const char *name[]) const =0
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
virtual void FillBin(Long64_t bin, Double_t w)=0
virtual Long64_t GetBin(const Double_t *x) const =0
void UpdateXStat(const Double_t *x, Double_t w=1.)
Definition THnBase.h:77
virtual Double_t GetBinContent(Long64_t bin, Int_t *idx=nullptr) const =0
void Multiply(const THnBase *h)
Multiply this histogram by histogram h this = this * h Note that if h has Sumw2 set,...
Definition THnBase.cxx:868
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:959
void Print(Option_t *option="") const override
Print a THnBase.
Definition THnBase.cxx:1446
Long64_t Fill(const char *name[], Double_t w=1.)
Definition THnBase.h:155
Double_t GetSumwx2(Int_t dim) const
Definition THnBase.h:214
Double_t GetWeightSum() const
Definition THnBase.h:139
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:554
void GetRandom(Double_t *rand, Bool_t subBinRandom=kTRUE)
Generate an n-dimensional random tuple based on the histogrammed distribution.
Definition THnBase.cxx:493
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:177
Long64_t Merge(TCollection *list)
Merge this with a list of THnBase's.
Definition THnBase.cxx:833
Double_t fTsumw
Total sum of weights.
Definition THnBase.h:49
virtual Long64_t GetBin(const Int_t *idx) const =0
Double_t ComputeIntegral()
Calculate the integral of the histogram.
Definition THnBase.cxx:1277
@ kInvalidInt
Definition THnBase.h:57
@ kNoInt
Definition THnBase.h:55
@ kValidInt
Definition THnBase.h:56
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:823
void SetBinError(Long64_t bin, Double_t e)
Definition THnBase.h:190
void CalculateErrors(Bool_t calc=kTRUE)
Calculate errors (or not if "calc" == kFALSE)
Definition THnBase.h:144
Int_t fNdimensions
Number of dimensions.
Definition THnBase.h:45
TH2D * Projection(Int_t yDim, Int_t xDim, Option_t *option="") const
Project all bins into a 2-dimensional histogram, keeping only axes "xDim" and "yDim".
Definition THnBase.h:237
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:193
void FillBinBase(Double_t w)
Increment the statistics due to filled weight "w",.
Definition THnBase.h:89
THnBase & operator=(const THnBase &other)
Definition THnBase.cxx:99
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
Bool_t HaveSkippedBin() const
Definition THnBase.h:345
THnIter(const THnBase *hist, Bool_t respectAxisRange=kFALSE)
Definition THnBase.h:332
Int_t GetCoord(Int_t dim) const
Definition THnBase.h:344
ROOT::Internal::THnBaseBinIter * fIter
Definition THnBase.h:349
~THnIter() override
Definition THnBase.cxx:1528
Bool_t RespectsAxisRange() const
Definition THnBase.h:346
Long64_t Next(Int_t *coord=nullptr)
Return the next bin's index.
Definition THnBase.h:340
A doubly linked list.
Definition TList.h:38
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
An array of TObjects.
Definition TObjArray.h:31
Mother of all ROOT objects.
Definition TObject.h:41
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:970
return c1
Definition legend1.C:41
Double_t x[n]
Definition legend1.C:17
TH1F * h1
Definition legend1.C:5
TF1 * f1
Definition legend1.C:11
return c2
Definition legend2.C:14
This file contains a specialised ROOT message handler to test for diagnostic in unit tests.
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:660