GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
GCube.h
Go to the documentation of this file.
1#ifndef GCUBE_H
2#define GCUBE_H
3
4#include "TH1.h"
5#include "TH2.h"
6#include "TH3.h"
7#include "TArrayF.h"
8#include "TArrayD.h"
9#include "TProfile.h"
10#include "TF1.h"
11#include "TRandom.h"
12
13class GCube : public TH1 {
14public:
15 GCube() = default;
16 GCube(const char* name, const char* title, Int_t nbins, Double_t low, Double_t up);
17 GCube(const char* name, const char* title, Int_t nbins, const Double_t* bins);
18 GCube(const char* name, const char* title, Int_t nbins, const Float_t* bins);
19 GCube(const GCube&);
20 GCube(GCube&&) noexcept;
21 GCube& operator=(const GCube&);
22 GCube& operator=(GCube&&) noexcept;
23
24 ~GCube() = default;
25
26 Int_t BufferEmpty(Int_t action = 0) override;
27 Int_t BufferFill(Double_t, Double_t) override { return -2; } // MayNotUse
28 virtual Int_t BufferFill(Double_t x, Double_t y, Double_t z, Double_t w);
29 void Copy(TObject& obj) const override;
30 Int_t Fill(Double_t) override; // MayNotUse
31 Int_t Fill(Double_t, Double_t) override { return Fill(0.); } // MayNotUse
32 Int_t Fill(const char*, Double_t) override { return Fill(0); } // MayNotUse
33 virtual Int_t Fill(Double_t, const char*, Double_t) { return Fill(0); } // MayNotUse
34 virtual Int_t Fill(const char*, Double_t, Double_t) { return Fill(0); } // MayNotUse
35 virtual Int_t Fill(const char*, const char*, Double_t) { return Fill(0); } // MayNotUse
36 virtual Int_t Fill(Double_t x, Double_t y, Double_t z);
37 virtual Int_t Fill(Double_t x, Double_t y, Double_t z, Double_t w);
38 virtual Int_t Fill(const char* namex, const char* namey, const char* namez, Double_t w);
39#if ROOT_VERSION_CODE < ROOT_VERSION(6, 24, 0)
40 void FillRandom(const char* fname, Int_t ntimes = 5000) override
41 {
42 FillRandom(fname, ntimes, nullptr);
43 }
44 void FillRandom(TH1* h, Int_t ntimes = 5000) override { FillRandom(h, ntimes, nullptr); }
45 void FillRandom(const char* fname, Int_t ntimes = 5000, TRandom* rng = nullptr);
46 void FillRandom(TH1* h, Int_t ntimes = 5000, TRandom* rng = nullptr);
47#elif ROOT_VERSION_CODE < ROOT_VERSION(6, 36, 0)
48 void FillRandom(const char* fname, Int_t ntimes = 5000, TRandom* rng = nullptr) override;
49 void FillRandom(TH1* h, Int_t ntimes = 5000, TRandom* rng = nullptr) override;
50#else
51 void FillRandom(const char* fname, Int_t ntimes = 5000, TRandom* rng = nullptr);
52 void FillRandom(TH1* h, Int_t ntimes = 5000, TRandom* rng = nullptr) override;
53#endif
54#if ROOT_VERSION_CODE < ROOT_VERSION(6, 18, 0)
55 Int_t FindFirstBinAbove(Double_t threshold = 0, Int_t axis = 1) const override
56 {
57 return FindFirstBinAbove(threshold, axis, 1, -1);
58 }
59 Int_t FindLastBinAbove(Double_t threshold = 0, Int_t axis = 1) const override { return FindLastBinAbove(threshold, axis, 1, -1); }
60 Int_t FindFirstBinAbove(Double_t threshold = 0, Int_t axis = 1, Int_t firstBin = 1, Int_t lastBin = -1) const;
61 Int_t FindLastBinAbove(Double_t threshold = 0, Int_t axis = 1, Int_t firstBin = 1, Int_t lastBin = -1) const;
62#else
63 Int_t FindFirstBinAbove(Double_t threshold = 0, Int_t axis = 1, Int_t firstBin = 1, Int_t lastBin = -1) const override;
64 Int_t FindLastBinAbove(Double_t threshold = 0, Int_t axis = 1, Int_t firstBin = 1, Int_t lastBin = -1) const override;
65#endif
66 virtual void FitSlicesZ(TF1* f1 = nullptr, Int_t binminx = 0, Int_t binmaxx = -1, Int_t binminy = 0,
67 Int_t binmaxy = -1, Int_t cut = 0, Option_t* option = "QNR");
68 Int_t GetBin(Int_t binx, Int_t biny = 0, Int_t binz = 0) const override;
69 virtual Double_t GetBinWithContent2(Double_t c, Int_t& binx, Int_t& biny, Int_t& binz, Int_t firstxbin = 1,
70 Int_t lastxbin = -1, Int_t firstybin = 1, Int_t lastybin = -1,
71 Int_t firstzbin = 1, Int_t lastzbin = -1, Double_t maxdiff = 0) const;
72 virtual Double_t GetCorrelationFactor(Int_t axis1 = 1, Int_t axis2 = 2) const;
73 virtual Double_t GetCovariance(Int_t axis1 = 1, Int_t axis2 = 2) const;
74 virtual void GetRandom3(Double_t& x, Double_t& y, Double_t& z);
75 void GetStats(Double_t* stats) const override;
76 Double_t Integral(Option_t* option = "") const override;
77 using TH1::Integral;
78 virtual Double_t Integral(Int_t firstxbin, Int_t lastxbin, Int_t firstybin, Int_t lastybin, Int_t firstzbin,
79 Int_t lastzbin, Option_t* option = "") const;
80 using TH1::IntegralAndError;
81 virtual Double_t IntegralAndError(Int_t firstxbin, Int_t lastxbin, Int_t firstybin, Int_t lastybin, Int_t firstzbin,
82 Int_t lastzbin, Double_t& error, Option_t* option = "") const;
83#if ROOT_VERSION_CODE < ROOT_VERSION(6, 20, 0)
84 Double_t Interpolate(Double_t) override;
85 Double_t Interpolate(Double_t, Double_t) override;
86 Double_t Interpolate(Double_t, Double_t, Double_t) override;
87#else
88 Double_t Interpolate(Double_t) const override;
89 Double_t Interpolate(Double_t, Double_t) const override;
90 Double_t Interpolate(Double_t, Double_t, Double_t) const override;
91#endif
92 Double_t KolmogorovTest(const TH1* h2, Option_t* option = "") const override;
93 Long64_t Merge(TCollection* list) override;
94 virtual TH1D* Projection(const char* name = "_pr", Int_t firstBiny = 0, Int_t lastBiny = -1, Int_t firstBinz = 0,
95 Int_t lastBinz = -1, Option_t* option = "") const;
96 void PutStats(Double_t* stats) override;
97 virtual GCube* Rebin3D(Int_t ngroup = 2, const char* newname = "");
98 void Reset(Option_t* option = "") override;
99 virtual void SetShowProjection(const char* option = "xy", Int_t nbins = 1); // *MENU*
100 TH1* ShowBackground(Int_t niter = 20, Option_t* option = "same") override;
101 Int_t ShowPeaks(Double_t sigma = 2, Option_t* option = "", Double_t threshold = 0.05) override; // *MENU*
102 void Smooth(Int_t ntimes = 1, Option_t* option = "") override; // *MENU*
103
104protected:
105 using TH1::DoIntegral;
106 Double_t DoIntegral(Int_t binx1, Int_t binx2, Int_t biny1, Int_t biny2, Int_t binz1, Int_t binz2, Double_t& error,
107 Option_t* option, Bool_t doError = kFALSE) const override;
108
109 void Matrix(TH2* val) { fMatrix = val; }
110 TH2* Matrix() { return fMatrix; }
111
112private:
113 Double_t fTsumwy{0}; // Total Sum of weight*Y
114 Double_t fTsumwy2{0}; // Total Sum of weight*Y*Y
115 Double_t fTsumwxy{0}; // Total Sum of weight*X*Y
116 Double_t fTsumwz{0}; // Total Sum of weight*Z
117 Double_t fTsumwz2{0}; // Total Sum of weight*Z*Z
118 Double_t fTsumwxz{0}; // Total Sum of weight*X*Z
119 Double_t fTsumwyz{0}; // Total Sum of weight*Y*Z
120 TH2* fMatrix{nullptr}; //!<! Transient pointer to the 2D-Matrix used in Draw() or GetMatrix()
121
122 /// /cond CLASSIMP
123 ClassDefOverride(GCube, 1) // NOLINT(readability-else-after-return)
124 /// /endcond
125};
126
127class GCubeF : public GCube, public TArrayF {
128public:
129 GCubeF();
130 GCubeF(const char* name, const char* title, Int_t nbins, Double_t low, Double_t up);
131 GCubeF(const char* name, const char* title, Int_t nbins, const Double_t* bins);
132 GCubeF(const char* name, const char* title, Int_t nbins, const Float_t* bins);
133 GCubeF(const GCubeF&);
134 GCubeF(GCubeF&&) noexcept;
136
137 TH2F* GetMatrix(bool force = false);
138
139 void AddBinContent(Int_t bin) override { ++fArray[bin]; }
140 void AddBinContent(Int_t bin, Double_t w) override { fArray[bin] += static_cast<Float_t>(w); }
141 void Copy(TObject& rh) const override;
142 void Draw(Option_t* option = "") override { GetMatrix()->Draw(option); }
143 TH1* DrawCopy(Option_t* option = "", const char* name_postfix = "_copy") const override;
144 Double_t GetBinContent(Int_t bin) const override;
145 Double_t GetBinContent(Int_t bin, Int_t) const override { return GetBinContent(bin); }
146 Double_t GetBinContent(Int_t binx, Int_t biny, Int_t binz) const override
147 {
148 return GetBinContent(GetBin(binx, biny, binz));
149 }
150 void Reset(Option_t* option = "") override;
151 Double_t RetrieveBinContent(Int_t bin) const override { return static_cast<Double_t>(fArray[bin]); }
152 void SetBinContent(Int_t bin, Double_t content) override;
153 void SetBinContent(Int_t bin, Int_t, Double_t content) override { SetBinContent(bin, content); }
154 void SetBinContent(Int_t binx, Int_t biny, Int_t binz, Double_t content) override
155 {
156 SetBinContent(GetBin(binx, biny, binz), content);
157 }
158 void SetBinsLength(Int_t n = -1) override;
159 void UpdateBinContent(Int_t bin, Double_t content) override { fArray[bin] = static_cast<Float_t>(content); }
160 GCubeF& operator=(const GCubeF& h1);
161 GCubeF& operator=(GCubeF&&) noexcept;
162 friend GCubeF operator*(Float_t c1, GCubeF& h1);
163 friend GCubeF operator*(GCubeF& h1, Float_t c1) { return operator*(c1, h1); }
164 friend GCubeF operator+(GCubeF& h1, GCubeF& h2);
165 friend GCubeF operator-(GCubeF& h1, GCubeF& h2);
166 friend GCubeF operator*(GCubeF& h1, GCubeF& h2);
167 friend GCubeF operator/(GCubeF& h1, GCubeF& h2);
168
169 /// /cond CLASSIMP
170 ClassDefOverride(GCubeF, 1) // NOLINT(readability-else-after-return)
171 /// /endcond
172};
173
174class GCubeD : public GCube, public TArrayD {
175public:
176 GCubeD();
177 GCubeD(const char* name, const char* title, Int_t nbins, Double_t low, Double_t up);
178 GCubeD(const char* name, const char* title, Int_t nbins, const Double_t* bins);
179 GCubeD(const char* name, const char* title, Int_t nbins, const Float_t* bins);
180 GCubeD(const GCubeD&);
181 GCubeD(GCubeD&&) noexcept;
183
184 TH2D* GetMatrix(bool force = false);
185
186 void AddBinContent(Int_t bin) override { ++fArray[bin]; }
187 void AddBinContent(Int_t bin, Double_t w) override { fArray[bin] += w; }
188 void Copy(TObject& rh) const override;
189 TH1* DrawCopy(Option_t* option = "", const char* name_postfix = "_copy") const override;
190 void Draw(Option_t* option = "") override { GetMatrix()->Draw(option); }
191 Double_t GetBinContent(Int_t bin) const override;
192 Double_t GetBinContent(Int_t bin, Int_t) const override { return GetBinContent(bin); }
193 Double_t GetBinContent(Int_t binx, Int_t biny, Int_t binz) const override
194 {
195 return GetBinContent(GetBin(binx, biny, binz));
196 }
197 void Reset(Option_t* option = "") override;
198 Double_t RetrieveBinContent(Int_t bin) const override { return static_cast<Double_t>(fArray[bin]); }
199 void SetBinContent(Int_t bin, Double_t content) override;
200 void SetBinContent(Int_t bin, Int_t, Double_t content) override { SetBinContent(bin, content); }
201 void SetBinContent(Int_t binx, Int_t biny, Int_t binz, Double_t content) override
202 {
203 SetBinContent(GetBin(binx, biny, binz), content);
204 }
205 void SetBinsLength(Int_t n = -1) override;
206 void UpdateBinContent(Int_t bin, Double_t content) override { fArray[bin] = content; }
207 GCubeD& operator=(const GCubeD& h1);
208 GCubeD& operator=(GCubeD&& h1) noexcept;
209 friend GCubeD operator*(Float_t c1, GCubeD& h1);
210 friend GCubeD operator*(GCubeD& h1, Float_t c1) { return operator*(c1, h1); }
211 friend GCubeD operator+(GCubeD& h1, GCubeD& h2);
212 friend GCubeD operator-(GCubeD& h1, GCubeD& h2);
213 friend GCubeD operator*(GCubeD& h1, GCubeD& h2);
214 friend GCubeD operator/(GCubeD& h1, GCubeD& h2);
215
216 /// /cond CLASSIMP
217 ClassDefOverride(GCubeD, 1) // NOLINT(readability-else-after-return)
218 /// /endcond
219};
220#endif
Double_t GetBinContent(Int_t binx, Int_t biny, Int_t binz) const override
Definition GCube.h:193
friend GCubeD operator*(GCubeD &h1, Float_t c1)
Definition GCube.h:210
TH1 * DrawCopy(Option_t *option="", const char *name_postfix="_copy") const override
Definition GCube.cxx:2930
void Reset(Option_t *option="") override
Definition GCube.cxx:2965
friend GCubeD operator+(GCubeD &h1, GCubeD &h2)
Definition GCube.cxx:3030
void SetBinsLength(Int_t n=-1) override
Definition GCube.cxx:2988
GCubeD()
Definition GCube.cxx:2851
void Draw(Option_t *option="") override
Definition GCube.h:190
void UpdateBinContent(Int_t bin, Double_t content) override
Definition GCube.h:206
GCubeD & operator=(const GCubeD &h1)
Definition GCube.cxx:3000
void AddBinContent(Int_t bin) override
Definition GCube.h:186
Double_t GetBinContent(Int_t bin, Int_t) const override
Definition GCube.h:192
TH2D * GetMatrix(bool force=false)
Definition GCube.cxx:2902
void Copy(TObject &rh) const override
Definition GCube.cxx:2925
void SetBinContent(Int_t bin, Double_t content) override
Definition GCube.cxx:2974
Double_t GetBinContent(Int_t bin) const override
Definition GCube.cxx:2947
friend GCubeD operator-(GCubeD &h1, GCubeD &h2)
Definition GCube.cxx:3040
friend GCubeD operator/(GCubeD &h1, GCubeD &h2)
Definition GCube.cxx:3060
friend GCubeD operator*(Float_t c1, GCubeD &h1)
Definition GCube.cxx:3020
void AddBinContent(Int_t bin, Double_t w) override
Definition GCube.h:187
void SetBinContent(Int_t bin, Int_t, Double_t content) override
Definition GCube.h:200
void SetBinContent(Int_t binx, Int_t biny, Int_t binz, Double_t content) override
Definition GCube.h:201
Double_t RetrieveBinContent(Int_t bin) const override
Definition GCube.h:198
GCubeF & operator=(const GCubeF &h1)
Definition GCube.cxx:2777
void SetBinContent(Int_t bin, Int_t, Double_t content) override
Definition GCube.h:153
void Draw(Option_t *option="") override
Definition GCube.h:142
void Reset(Option_t *option="") override
Definition GCube.cxx:2742
void AddBinContent(Int_t bin) override
Definition GCube.h:139
Double_t GetBinContent(Int_t binx, Int_t biny, Int_t binz) const override
Definition GCube.h:146
GCubeF()
Definition GCube.cxx:2629
friend GCubeF operator-(GCubeF &h1, GCubeF &h2)
Definition GCube.cxx:2817
TH1 * DrawCopy(Option_t *option="", const char *name_postfix="_copy") const override
Definition GCube.cxx:2706
friend GCubeF operator/(GCubeF &h1, GCubeF &h2)
Definition GCube.cxx:2837
Double_t RetrieveBinContent(Int_t bin) const override
Definition GCube.h:151
void SetBinsLength(Int_t n=-1) override
Definition GCube.cxx:2765
TH2F * GetMatrix(bool force=false)
Definition GCube.cxx:2680
Double_t GetBinContent(Int_t bin, Int_t) const override
Definition GCube.h:145
void UpdateBinContent(Int_t bin, Double_t content) override
Definition GCube.h:159
void SetBinContent(Int_t binx, Int_t biny, Int_t binz, Double_t content) override
Definition GCube.h:154
Double_t GetBinContent(Int_t bin) const override
Definition GCube.cxx:2723
friend GCubeF operator*(GCubeF &h1, Float_t c1)
Definition GCube.h:163
void Copy(TObject &rh) const override
Definition GCube.cxx:2701
friend GCubeF operator+(GCubeF &h1, GCubeF &h2)
Definition GCube.cxx:2807
void SetBinContent(Int_t bin, Double_t content) override
Definition GCube.cxx:2751
void AddBinContent(Int_t bin, Double_t w) override
Definition GCube.h:140
friend GCubeF operator*(Float_t c1, GCubeF &h1)
Definition GCube.cxx:2797
Definition GCube.h:13
TH2 * fMatrix
! Transient pointer to the 2D-Matrix used in Draw() or GetMatrix()
Definition GCube.h:120
Int_t Fill(Double_t) override
Definition GCube.cxx:309
virtual Double_t GetCovariance(Int_t axis1=1, Int_t axis2=2) const
Definition GCube.cxx:1042
Double_t fTsumwxz
Definition GCube.h:118
Double_t fTsumwz
Definition GCube.h:116
Int_t Fill(const char *, Double_t) override
Definition GCube.h:32
virtual GCube * Rebin3D(Int_t ngroup=2, const char *newname="")
Definition GCube.cxx:2097
Int_t Fill(Double_t, Double_t) override
Definition GCube.h:31
Int_t ShowPeaks(Double_t sigma=2, Option_t *option="", Double_t threshold=0.05) override
Definition GCube.cxx:2497
TH2 * Matrix()
Definition GCube.h:110
virtual void FitSlicesZ(TF1 *f1=nullptr, Int_t binminx=0, Int_t binmaxx=-1, Int_t binminy=0, Int_t binmaxy=-1, Int_t cut=0, Option_t *option="QNR")
Definition GCube.cxx:812
Double_t fTsumwy2
Definition GCube.h:114
virtual Int_t Fill(const char *, Double_t, Double_t)
Definition GCube.h:34
void Smooth(Int_t ntimes=1, Option_t *option="") override
Definition GCube.cxx:2510
Double_t Interpolate(Double_t) const override
Definition GCube.cxx:1286
virtual Double_t IntegralAndError(Int_t firstxbin, Int_t lastxbin, Int_t firstybin, Int_t lastybin, Int_t firstzbin, Int_t lastzbin, Double_t &error, Option_t *option="") const
Definition GCube.cxx:1270
Int_t GetBin(Int_t binx, Int_t biny=0, Int_t binz=0) const override
Definition GCube.cxx:909
virtual void SetShowProjection(const char *option="xy", Int_t nbins=1)
Definition GCube.cxx:2459
GCube()=default
Double_t fTsumwy
Definition GCube.h:113
virtual void GetRandom3(Double_t &x, Double_t &y, Double_t &z)
Definition GCube.cxx:1089
TH1 * ShowBackground(Int_t niter=20, Option_t *option="same") override
Definition GCube.cxx:2487
virtual Int_t Fill(const char *, const char *, Double_t)
Definition GCube.h:35
Int_t FindLastBinAbove(Double_t threshold=0, Int_t axis=1, Int_t firstBin=1, Int_t lastBin=-1) const override
Definition GCube.cxx:731
void PutStats(Double_t *stats) override
Definition GCube.cxx:2084
virtual Int_t Fill(Double_t, const char *, Double_t)
Definition GCube.h:33
Double_t Integral(Option_t *option="") const override
Definition GCube.cxx:1247
virtual Double_t GetCorrelationFactor(Int_t axis1=1, Int_t axis2=2) const
Definition GCube.cxx:1020
void GetStats(Double_t *stats) const override
Definition GCube.cxx:1138
void Reset(Option_t *option="") override
Definition GCube.cxx:2435
void Matrix(TH2 *val)
Definition GCube.h:109
void FillRandom(const char *fname, Int_t ntimes=5000, TRandom *rng=nullptr)
Definition GCube.cxx:558
Int_t BufferEmpty(Int_t action=0) override
Definition GCube.cxx:70
virtual TH1D * Projection(const char *name="_pr", Int_t firstBiny=0, Int_t lastBiny=-1, Int_t firstBinz=0, Int_t lastBinz=-1, Option_t *option="") const
Definition GCube.cxx:1832
Double_t DoIntegral(Int_t binx1, Int_t binx2, Int_t biny1, Int_t biny2, Int_t binz1, Int_t binz2, Double_t &error, Option_t *option, Bool_t doError=kFALSE) const override
Definition GCube.cxx:223
Long64_t Merge(TCollection *list) override
Definition GCube.cxx:1581
Double_t KolmogorovTest(const TH1 *h2, Option_t *option="") const override
Definition GCube.cxx:1369
Double_t fTsumwz2
Definition GCube.h:117
Int_t FindFirstBinAbove(Double_t threshold=0, Int_t axis=1, Int_t firstBin=1, Int_t lastBin=-1) const override
Definition GCube.cxx:682
Int_t BufferFill(Double_t, Double_t) override
Definition GCube.h:27
Double_t fTsumwyz
Definition GCube.h:119
void Copy(TObject &obj) const override
Definition GCube.cxx:212
Double_t fTsumwxy
Definition GCube.h:115
virtual Double_t GetBinWithContent2(Double_t c, Int_t &binx, Int_t &biny, Int_t &binz, Int_t firstxbin=1, Int_t lastxbin=-1, Int_t firstybin=1, Int_t lastybin=-1, Int_t firstzbin=1, Int_t lastzbin=-1, Double_t maxdiff=0) const
Definition GCube.cxx:946