GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
GHSym.cxx
Go to the documentation of this file.
1#include "GHSym.h"
2
3#include "TROOT.h"
4#include "THashList.h"
5#include "THLimitsFinder.h"
6#include "TVirtualHistPainter.h"
7#include "TObjString.h"
8#include "TVirtualPad.h"
9#include "TMath.h"
10#include "TPad.h"
11#include "TF1.h"
12#include "TF2.h"
13#include "TRandom.h"
14#include "TClass.h"
15
16#include <iostream>
17
18// Internal exceptions for the CheckConsistency method
19class DifferentDimension : public std::exception {
20};
21class DifferentNumberOfBins : public std::exception {
22};
23class DifferentAxisLimits : public std::exception {
24};
25class DifferentBinLimits : public std::exception {
26};
27class DifferentLabels : public std::exception {
28};
29
31{
32 fDimension = 2;
33}
34
35// we have to repeat the code from the default constructor here, because calling the TH1 constructor and the GHSym
36// default constructor gives an error
37GHSym::GHSym(const char* name, const char* title, Int_t nbins, Double_t low, Double_t up)
38 : TH1(name, title, nbins, low, up)
39{
40 fDimension = 2;
41 fYaxis.Set(nbins, low, up);
42 // TH1 constructor sets fNcells to nbins+2
43 // we need (nbins+2)*((nbins+2)+1)/2 cells
44 fNcells = (fNcells * (nbins + 3)) / 2;
45}
46
47GHSym::GHSym(const char* name, const char* title, Int_t nbins, const Double_t* bins)
48 : TH1(name, title, nbins, bins)
49{
50 fDimension = 2;
51 fYaxis.Set(nbins, bins);
52 // TH1 constructor sets fNcells to nbins+2
53 // we need (nbins+2)*((nbins+2)+1)/2 cells
54 fNcells = (fNcells * (nbins + 3)) / 2;
55}
56
57GHSym::GHSym(const char* name, const char* title, Int_t nbins, const Float_t* bins)
58 : TH1(name, title, nbins, bins)
59{
60 fDimension = 2;
61 fYaxis.Set(nbins, bins);
62 // TH1 constructor sets fNcells to nbins+2
63 // we need (nbins+2)*((nbins+2)+1)/2 cells
64 fNcells = (fNcells * (nbins + 3)) / 2;
65}
66
67GHSym::GHSym(const GHSym& rhs) : TH1() // NOLINT(readability-redundant-member-init)
68{
69 rhs.Copy(*this);
70}
71
72GHSym::GHSym(GHSym&& rhs) noexcept : TH1() // NOLINT(readability-redundant-member-init)
73{
74 rhs.Copy(*this);
75}
76
77Int_t GHSym::BufferEmpty(Int_t action)
78{
79 /// Fill histogram with all entries in the buffer.
80 /// action = -1 histogram is reset and refilled from the buffer (called by THistPainter::Paint)
81 /// action = 0 histogram is filled from the buffer
82 /// action = 1 histogram is filled and buffer is deleted
83 /// The buffer is automatically deleted when the number of entries
84 /// in the buffer is greater than the number of entries in the histogram
85 if(fBuffer == nullptr) {
86 return 0;
87 }
88
89 auto nbEntries = static_cast<Int_t>(fBuffer[0]);
90 if(nbEntries == 0) {
91 return 0;
92 }
93 if(nbEntries < 0 && action == 0) {
94 return 0; // histogram has been already filled from the buffer
95 }
96 Double_t* buffer = fBuffer;
97 if(nbEntries < 0) {
98 nbEntries = -nbEntries;
99 fBuffer = nullptr;
100 Reset("ICES");
101 fBuffer = buffer;
102 }
103
104 if(CanExtendAllAxes() || fXaxis.GetXmax() <= fXaxis.GetXmin() || fYaxis.GetXmax() <= fYaxis.GetXmin()) {
105 // find min, max of entries in buffer
106 // for the symmetric matrix x- and y-range are the same
107 Double_t min = fBuffer[2];
108 Double_t max = min;
109 if(fBuffer[3] < min) {
110 min = fBuffer[3];
111 }
112 if(fBuffer[3] > max) {
113 max = fBuffer[3];
114 }
115 for(Int_t i = 1; i < nbEntries; ++i) {
116 Double_t x = fBuffer[3 * i + 2];
117 if(x < min) {
118 min = x;
119 }
120 if(x > max) {
121 max = x;
122 }
123 Double_t y = fBuffer[3 * i + 3];
124 if(y < min) {
125 min = y;
126 }
127 if(y > max) {
128 max = y;
129 }
130 }
131 if(fXaxis.GetXmax() <= fXaxis.GetXmin() || fYaxis.GetXmax() <= fYaxis.GetXmin()) {
132 THLimitsFinder::GetLimitsFinder()->FindGoodLimits(this, min, max, min, max);
133 } else {
134 fBuffer = nullptr;
135 Int_t keep = fBufferSize;
136 fBufferSize = 0;
137 if(min < fXaxis.GetXmin()) {
138 RebinAxis(min, &fXaxis);
139 }
140 if(max >= fXaxis.GetXmax()) {
141 RebinAxis(max, &fXaxis);
142 }
143 if(min < fYaxis.GetXmin()) {
144 RebinAxis(min, &fYaxis);
145 }
146 if(max >= fYaxis.GetXmax()) {
147 RebinAxis(max, &fYaxis);
148 }
149 fBuffer = buffer;
150 fBufferSize = keep;
151 }
152 }
153
154 fBuffer = nullptr;
155 for(Int_t i = 0; i < nbEntries; ++i) {
156 Fill(buffer[3 * i + 2], buffer[3 * i + 3], buffer[3 * i + 1]);
157 }
158 fBuffer = buffer;
159
160 if(action > 0) {
161 delete[] fBuffer;
162 fBuffer = nullptr;
163 fBufferSize = 0;
164 } else {
165 if(nbEntries == static_cast<Int_t>(fEntries)) {
166 fBuffer[0] = -nbEntries;
167 } else {
168 fBuffer[0] = 0;
169 }
170 }
171 return nbEntries;
172}
173
174Int_t GHSym::BufferFill(Double_t x, Double_t y, Double_t w)
175{
176 if(fBuffer == nullptr) {
177 return -3;
178 }
179
180 auto nbEntries = static_cast<Int_t>(fBuffer[0]);
181 if(nbEntries < 0) {
182 nbEntries = -nbEntries;
183 fBuffer[0] = nbEntries;
184 if(fEntries > 0) {
185 Double_t* buffer = fBuffer;
186 fBuffer = nullptr;
187 Reset("ICES");
188 fBuffer = buffer;
189 }
190 }
191 if(3 * nbEntries + 3 >= fBufferSize) {
192 BufferEmpty(1);
193 return Fill(x, y, w);
194 }
195 fBuffer[3 * nbEntries + 1] = w;
196 fBuffer[3 * nbEntries + 2] = x;
197 fBuffer[3 * nbEntries + 3] = y;
198 fBuffer[0] += 1;
199
200 return -3;
201}
202
203void GHSym::Copy(TObject& obj) const
204{
205 // Copy.
206
207 TH1::Copy(obj);
208 static_cast<GHSym&>(obj).fTsumwy = fTsumwy;
209 static_cast<GHSym&>(obj).fTsumwy2 = fTsumwy2;
210 static_cast<GHSym&>(obj).fTsumwxy = fTsumwxy;
211 static_cast<GHSym&>(obj).fMatrix = nullptr;
212}
213
214Double_t GHSym::DoIntegral(Int_t binx1, Int_t binx2, Int_t biny1, Int_t biny2, Double_t& error, Option_t* option, Bool_t doError) const
215{
216 // internal function compute integral and optionally the error between the limits
217 // specified by the bin number values working for all histograms (1D, 2D and 3D)
218
219 Int_t nbinsx = GetNbinsX();
220 if(binx1 < 0) {
221 binx1 = 0;
222 }
223 if(binx2 > nbinsx + 1 || binx2 < binx1) {
224 binx2 = nbinsx + 1;
225 }
226 if(GetDimension() > 1) {
227 Int_t nbinsy = GetNbinsY();
228 if(biny1 < 0) {
229 biny1 = 0;
230 }
231 if(biny2 > nbinsy + 1 || biny2 < biny1) {
232 biny2 = nbinsy + 1;
233 }
234 } else {
235 biny1 = 0;
236 biny2 = 0;
237 }
238
239 // - Loop on bins in specified range
240 TString opt = option;
241 opt.ToLower();
242 Bool_t width = kFALSE;
243 if(opt.Contains("width")) {
244 width = kTRUE;
245 }
246
247 Double_t dx = 1.;
248 Double_t dy = 1.;
249 Double_t integral = 0;
250 Double_t igerr2 = 0;
251 for(Int_t binx = binx1; binx <= binx2; ++binx) {
252 if(width) {
253 dx = fXaxis.GetBinWidth(binx);
254 }
255 for(Int_t biny = biny1; biny <= biny2; ++biny) {
256 if(width) {
257 dy = fYaxis.GetBinWidth(biny);
258 }
259 Int_t bin = GetBin(binx, biny);
260 if(width) {
261 integral += GetBinContent(bin) * dx * dy;
262 } else {
263 integral += GetBinContent(bin);
264 }
265 if(doError) {
266 if(width) {
267 igerr2 += GetBinError(bin) * GetBinError(bin) * dx * dx * dy * dy;
268 } else {
269 igerr2 += GetBinError(bin) * GetBinError(bin);
270 }
271 }
272 }
273 }
274
275 if(doError) {
276 error = TMath::Sqrt(igerr2);
277 }
278 return integral;
279}
280
281Int_t GHSym::Fill(Double_t)
282{
283 // Invalid Fill method
284 Error("Fill", "Invalid signature - do nothing");
285 return -1;
286}
287
288Int_t GHSym::Fill(Double_t x, Double_t y)
289{
290 /// Increment cell defined by x,y by 1.
291 if(fBuffer != nullptr) {
292 return BufferFill(x, y, 1);
293 }
294
295 Int_t binx = 0;
296 Int_t biny = 0;
297 fEntries++;
298 if(y <= x) {
299 binx = fXaxis.FindBin(x);
300 biny = fYaxis.FindBin(y);
301 } else {
302 binx = fXaxis.FindBin(y);
303 biny = fYaxis.FindBin(x);
304 }
305 if(binx < 0 || biny < 0) {
306 return -1;
307 }
308 Int_t bin = biny * (2 * fXaxis.GetNbins() - biny + 3) / 2 + binx;
309 AddBinContent(bin);
310 if(fSumw2.fN != 0) {
311 ++fSumw2.fArray[bin];
312 }
313 if(binx == 0 || binx > fXaxis.GetNbins()) {
314 if(!fgStatOverflows) {
315 return -1;
316 }
317 }
318 if(biny == 0 || biny > fYaxis.GetNbins()) {
319 if(!fgStatOverflows) {
320 return -1;
321 }
322 }
323 // not sure if these summed weights are calculated correct
324 // as of now this is the method used in TH2
325 ++fTsumw;
326 ++fTsumw2;
327 fTsumwx += x;
328 fTsumwx2 += x * x;
329 fTsumwy += y;
330 fTsumwy2 += y * y;
331 fTsumwxy += x * y;
332 return bin;
333}
334
335Int_t GHSym::Fill(Double_t x, Double_t y, Double_t w)
336{
337 /// Increment cell defined by x,y by 1.
338 if(fBuffer != nullptr) {
339 return BufferFill(x, y, 1);
340 }
341
342 Int_t binx = 0;
343 Int_t biny = 0;
344 fEntries++;
345 if(y <= x) {
346 binx = fXaxis.FindBin(x);
347 biny = fYaxis.FindBin(y);
348 } else {
349 binx = fXaxis.FindBin(y);
350 biny = fYaxis.FindBin(x);
351 }
352 if(binx < 0 || biny < 0) {
353 return -1;
354 }
355 Int_t bin = biny * (2 * fXaxis.GetNbins() - biny + 3) / 2 + binx;
356 AddBinContent(bin, w);
357 if(fSumw2.fN != 0) {
358 fSumw2.fArray[bin] += w * w;
359 }
360 if(binx == 0 || binx > fXaxis.GetNbins()) {
361 if(!fgStatOverflows) {
362 return -1;
363 }
364 }
365 if(biny == 0 || biny > fYaxis.GetNbins()) {
366 if(!fgStatOverflows) {
367 return -1;
368 }
369 }
370 // not sure if these summed weights are calculated correct
371 // as of now this is the method used in TH2
372 fTsumw += w;
373 fTsumw2 += w * w;
374 fTsumwx += w * x;
375 fTsumwx2 += w * x * x;
376 fTsumwy += w * y;
377 fTsumwy2 += w * y * y;
378 fTsumwxy += w * x * y;
379 return bin;
380}
381
382Int_t GHSym::Fill(const char* namex, const char* namey, Double_t w)
383{
384 // Increment cell defined by namex,namey by a weight w
385 //
386 // if x or/and y is less than the low-edge of the corresponding axis first bin,
387 // the Underflow cell is incremented.
388 // if x or/and y is greater than the upper edge of corresponding axis last bin,
389 // the Overflow cell is incremented.
390 //
391 // If the storage of the sum of squares of weights has been triggered,
392 // via the function Sumw2, then the sum of the squares of weights is incremented
393 // by w^2 in the cell corresponding to x,y.
394 //
395
396 Int_t binx = 0;
397 Int_t biny = 0;
398 fEntries++;
399 binx = fXaxis.FindBin(namex);
400 biny = fYaxis.FindBin(namey);
401 if(binx < 0 || biny < 0) {
402 return -1;
403 }
404 if(biny >= binx) {
405 std::swap(binx, biny);
406 }
407 Int_t bin = biny * (2 * fXaxis.GetNbins() - biny + 3) / 2 + binx;
408 AddBinContent(bin, w);
409 if(fSumw2.fN != 0) {
410 fSumw2.fArray[bin] += w * w;
411 }
412 if(binx == 0 || binx > fXaxis.GetNbins()) {
413 return -1;
414 }
415 if(biny == 0 || biny > fYaxis.GetNbins()) {
416 return -1;
417 }
418 Double_t x = fXaxis.GetBinCenter(binx);
419 Double_t y = fYaxis.GetBinCenter(biny);
420 fTsumw += w;
421 fTsumw2 += w * w;
422 fTsumwx += w * x;
423 fTsumwx2 += w * x * x;
424 fTsumwy += w * y;
425 fTsumwy2 += w * y * y;
426 fTsumwxy += w * x * y;
427 return bin;
428}
429
430void GHSym::FillN(Int_t ntimes, const Double_t* x, const Double_t* y, const Double_t* w, Int_t stride)
431{
432 //*-*-*-*-*-*-*Fill a 2-D histogram with an array of values and weights*-*-*-*
433 //*-* ========================================================
434 //*-*
435 //*-* ntimes: number of entries in arrays x and w (array size must be ntimes*stride)
436 //*-* x: array of x values to be histogrammed
437 //*-* y: array of y values to be histogrammed
438 //*-* w: array of weights
439 //*-* stride: step size through arrays x, y and w
440 //*-*
441 //*-* If the storage of the sum of squares of weights has been triggered,
442 //*-* via the function Sumw2, then the sum of the squares of weights is incremented
443 //*-* by w[i]^2 in the cell corresponding to x[i],y[i].
444 //*-* if w is NULL each entry is assumed a weight=1
445 //*-*
446 //*-* NB: function only valid for a TH2x object
447 //*-*
448 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
449 ntimes *= stride;
450 for(int i = 0; i < ntimes; i += stride) {
451 if(w != nullptr) {
452 Fill(x[i], y[i], w[i]);
453 } else {
454 Fill(x[i], y[i]);
455 }
456 }
457}
458
459void GHSym::FillRandom(const char* fname, Int_t ntimes, TRandom* rng)
460{
461 //*-*-*-*-*-*-*Fill histogram following distribution in function fname*-*-*-*
462 //*-* =======================================================
463 //*-*
464 //*-* The distribution contained in the function fname (TF2) is integrated
465 //*-* over the channel contents.
466 //*-* It is normalized to 1.
467 //*-* Getting one random number implies:
468 //*-* - Generating a random number between 0 and 1 (say r1)
469 //*-* - Look in which bin in the normalized integral r1 corresponds to
470 //*-* - Fill histogram channel
471 //*-* ntimes random numbers are generated
472 //*-*
473 //*-* One can also call TF2::GetRandom2 to get a random variate from a function.
474 //*-*
475 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-**-*-*-*-*-*-*-*
476
477 //*-*- Search for fname in the list of ROOT defined functions
478 TObject* fobj = gROOT->GetFunction(fname);
479 if(fobj == nullptr) {
480 Error("FillRandom", "Unknown function: %s", fname);
481 return;
482 }
483 TF2* f1 = static_cast<TF2*>(fobj);
484 if(f1 == nullptr) {
485 Error("FillRandom", "Function: %s is not a TF2", fname);
486 return;
487 }
488
489 //*-*- Allocate temporary space to store the integral and compute integral
490 Int_t nbinsx = GetNbinsX();
491 Int_t nbinsy = GetNbinsY();
492 Int_t nbins = nbinsx * nbinsy;
493
494 auto* integral = new Double_t[nbins + 1];
495 Int_t ibin = 0;
496 integral[ibin] = 0;
497 for(Int_t biny = 1; biny <= nbinsy; ++biny) {
498 for(Int_t binx = 1; binx <= nbinsx; ++binx) {
499 ++ibin;
500 Double_t fint = f1->Integral(fXaxis.GetBinLowEdge(binx), fXaxis.GetBinUpEdge(binx), fYaxis.GetBinLowEdge(biny),
501 fYaxis.GetBinUpEdge(biny));
502 integral[ibin] = integral[ibin - 1] + fint;
503 }
504 }
505
506 //*-*- Normalize integral to 1
507 if(integral[nbins] == 0) {
508 delete[] integral;
509 Error("FillRandom", "Integral = zero");
510 return;
511 }
512 for(Int_t bin = 1; bin <= nbins; ++bin) {
513 integral[bin] /= integral[nbins];
514 }
515
516 //*-*--------------Start main loop ntimes
517 for(int loop = 0; loop < ntimes; ++loop) {
518 Double_t r1 = (rng != nullptr) ? rng->Rndm(loop) : gRandom->Rndm(loop);
519 ibin = static_cast<int>(TMath::BinarySearch(nbins, &integral[0], r1));
520 Int_t biny = ibin / nbinsx;
521 Int_t binx = 1 + ibin - nbinsx * biny;
522 ++biny;
523 Fill(fXaxis.GetBinCenter(binx), fYaxis.GetBinCenter(biny));
524 }
525 delete[] integral;
526}
527
528#if ROOT_VERSION_CODE < ROOT_VERSION(6, 24, 0)
529void GHSym::FillRandom(TH1* h, Int_t ntimes, TRandom*)
530#else
531void GHSym::FillRandom(TH1* hist, Int_t ntimes, TRandom* rng)
532#endif
533{
534 //*-*-*-*-*-*-*Fill histogram following distribution in histogram h*-*-*-*
535 //*-* ====================================================
536 //*-*
537 //*-* The distribution contained in the histogram h (TH2) is integrated
538 //*-* over the channel contents.
539 //*-* It is normalized to 1.
540 //*-* Getting one random number implies:
541 //*-* - Generating a random number between 0 and 1 (say r1)
542 //*-* - Look in which bin in the normalized integral r1 corresponds to
543 //*-* - Fill histogram channel
544 //*-* ntimes random numbers are generated
545 //*-*
546 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-**-*-*-*-*-*-*-*
547
548 if(hist == nullptr) {
549 Error("FillRandom", "Null histogram");
550 return;
551 }
552 if(fDimension != hist->GetDimension()) {
553 Error("FillRandom", "Histograms with different dimensions");
554 return;
555 }
556
557 if(hist->ComputeIntegral() == 0) {
558 return;
559 }
560
561 Double_t x = 0.;
562 Double_t y = 0.;
563 TH2* h2 = static_cast<TH2*>(hist);
564 for(int loop = 0; loop < ntimes; ++loop) {
565#if ROOT_VERSION_CODE < ROOT_VERSION(6, 24, 0)
566 h2->GetRandom2(x, y);
567#else
568 h2->GetRandom2(x, y, rng);
569#endif
570 Fill(x, y);
571 }
572}
573
574Int_t GHSym::FindFirstBinAbove(Double_t threshold, Int_t axis, Int_t firstBin, Int_t lastBin) const
575{
576 // find first bin with content > threshold for axis (1=x, 2=y, 3=z)
577 // if no bins with content > threshold is found the function returns -1.
578
579 if(axis < 1 || axis > 2) {
580 Warning("FindFirstBinAbove", "Invalid axis number : %d, axis x assumed\n", axis);
581 axis = 1;
582 }
583 Int_t nbinsx = fXaxis.GetNbins();
584 if(lastBin > firstBin && lastBin < nbinsx) { nbinsx = lastBin; }
585 Int_t nbinsy = fYaxis.GetNbins();
586 if(lastBin > firstBin && lastBin < nbinsy) { nbinsy = lastBin; }
587 if(axis == 1) {
588 for(Int_t binx = firstBin; binx <= nbinsx; ++binx) {
589 for(Int_t biny = firstBin; biny <= nbinsy; ++biny) {
590 if(GetBinContent(binx, biny) > threshold) {
591 return binx;
592 }
593 }
594 }
595 } else {
596 for(Int_t biny = firstBin; biny <= nbinsy; ++biny) {
597 for(Int_t binx = firstBin; binx <= nbinsx; ++binx) {
598 if(GetBinContent(binx, biny) > threshold) {
599 return biny;
600 }
601 }
602 }
603 }
604 return -1;
605}
606
607Int_t GHSym::FindLastBinAbove(Double_t threshold, Int_t axis, Int_t firstBin, Int_t lastBin) const
608{
609 // find last bin with content > threshold for axis (1=x, 2=y, 3=z)
610 // if no bins with content > threshold is found the function returns -1.
611
612 if(axis < 1 || axis > 2) {
613 Warning("FindLastBinAbove", "Invalid axis number : %d, axis x assumed\n", axis);
614 axis = 1;
615 }
616 Int_t nbinsx = fXaxis.GetNbins();
617 if(lastBin > firstBin && lastBin < nbinsx) { nbinsx = lastBin; }
618 Int_t nbinsy = fYaxis.GetNbins();
619 if(lastBin > firstBin && lastBin < nbinsy) { nbinsy = lastBin; }
620 if(axis == 1) {
621 for(Int_t binx = nbinsx; binx >= firstBin; --binx) {
622 for(Int_t biny = firstBin; biny <= nbinsy; ++biny) {
623 if(GetBinContent(binx, biny) > threshold) {
624 return binx;
625 }
626 }
627 }
628 } else {
629 for(Int_t biny = nbinsy; biny >= firstBin; --biny) {
630 for(Int_t binx = firstBin; binx <= nbinsx; ++binx) {
631 if(GetBinContent(binx, biny) > threshold) {
632 return biny;
633 }
634 }
635 }
636 }
637 return -1;
638}
639
640void GHSym::FitSlices(TF1* f1, Int_t firstbin, Int_t lastbin, Int_t cut, Option_t* option, TObjArray* arr)
641{
642 // Project slices along X in case of a 2-D histogram, then fit each slice
643 // with function f1 and make a histogram for each fit parameter
644 // Only bins along Y between firstbin and lastbin are considered.
645 // By default (firstbin == 0, lastbin == -1), all bins in y including
646 // over- and underflows are taken into account.
647 // If f1=0, a gaussian is assumed
648 // Before invoking this function, one can set a subrange to be fitted along X
649 // via f1->SetRange(xmin,xmax)
650 // The argument option (default="QNR") can be used to change the fit options.
651 // "Q" means Quiet mode
652 // "N" means do not show the result of the fit
653 // "R" means fit the function in the specified function range
654 // "G2" merge 2 consecutive bins along X
655 // "G3" merge 3 consecutive bins along X
656 // "G4" merge 4 consecutive bins along X
657 // "G5" merge 5 consecutive bins along X
658 // "S" sliding merge: merge n consecutive bins along X accordingly to what Gn is given.
659 // It makes sense when used together with a Gn option
660 //
661 // The generated histograms are returned by adding them to arr, if arr is not NULL.
662 // arr's SetOwner() is called, to signal that it is the user's respponsability to
663 // delete the histograms, possibly by deleting the arrary.
664 // TObjArray aSlices;
665 // h2->FitSlicesX(func, 0, -1, 0, "QNR", &aSlices);
666 // will already delete the histograms once aSlice goes out of scope. aSlices will
667 // contain the histogram for the i-th parameter of the fit function at aSlices[i];
668 // aSlices[n] (n being the number of parameters) contains the chi2 distribution of
669 // the fits.
670 //
671 // If arr is NULL, the generated histograms are added to the list of objects
672 // in the current directory. It is the user's responsability to delete
673 // these histograms.
674 //
675 // Example: Assume a 2-d histogram h2
676 // Root > h2->FitSlicesX(); produces 4 TH1D histograms
677 // with h2_0 containing parameter 0(Constant) for a Gaus fit
678 // of each bin in Y projected along X
679 // with h2_1 containing parameter 1(Mean) for a gaus fit
680 // with h2_2 containing parameter 2(RMS) for a gaus fit
681 // with h2_chi2 containing the chisquare/number of degrees of freedom for a gaus fit
682 //
683 // Root > h2->FitSlicesX(0,15,22,10);
684 // same as above, but only for bins 15 to 22 along Y
685 // and only for bins in Y for which the corresponding projection
686 // along X has more than cut bins filled.
687 //
688 // NOTE: To access the generated histograms in the current directory, do eg:
689 // TH1D *h2_1 = (TH1D*)gDirectory->Get("h2_1");
690
691 Int_t nbins = fYaxis.GetNbins();
692 if(firstbin < 0) {
693 firstbin = 0;
694 }
695 if(lastbin < 0 || lastbin > nbins + 1) {
696 lastbin = nbins + 1;
697 }
698 if(lastbin < firstbin) {
699 firstbin = 0;
700 lastbin = nbins + 1;
701 }
702 TString opt = option;
703 opt.ToLower();
704 Int_t ngroup = 1;
705 if(opt.Contains("g2")) {
706 ngroup = 2;
707 opt.ReplaceAll("g2", "");
708 }
709 if(opt.Contains("g3")) {
710 ngroup = 3;
711 opt.ReplaceAll("g3", "");
712 }
713 if(opt.Contains("g4")) {
714 ngroup = 4;
715 opt.ReplaceAll("g4", "");
716 }
717 if(opt.Contains("g5")) {
718 ngroup = 5;
719 opt.ReplaceAll("g5", "");
720 }
721
722 // implement option S sliding merge for each bin using in conjunction with a given Gn
723 Int_t nstep = ngroup;
724 if(opt.Contains("s")) {
725 nstep = 1;
726 }
727
728 // default is to fit with a gaussian
729 if(f1 == nullptr) {
730 f1 = static_cast<TF1*>(gROOT->GetFunction("gaus"));
731 if(f1 == nullptr) {
732 f1 = new TF1("gaus", "gaus", fXaxis.GetXmin(), fXaxis.GetXmax());
733 } else {
734 f1->SetRange(fXaxis.GetXmin(), fXaxis.GetXmax());
735 }
736 }
737 Int_t npar = f1->GetNpar();
738 if(npar <= 0) {
739 return;
740 }
741 auto* parsave = new Double_t[npar];
742 f1->GetParameters(parsave);
743
744 if(arr != nullptr) {
745 arr->SetOwner();
746 arr->Expand(npar + 1);
747 }
748
749 // Create one histogram for each function parameter
750 auto** hlist = new TH1D*[npar];
751 auto* name = new char[2000];
752 auto* title = new char[2000];
753 const TArrayD* bins = fYaxis.GetXbins();
754 for(Int_t ipar = 0; ipar < npar; ++ipar) {
755 snprintf(name, 2000, "%s_%d", GetName(), ipar);
756 snprintf(title, 2000, "Fitted value of par[%d]=%s", ipar, f1->GetParName(ipar));
757 delete gDirectory->FindObject(name);
758 if(bins->fN == 0) {
759 hlist[ipar] = new TH1D(name, title, nbins, fYaxis.GetXmin(), fYaxis.GetXmax());
760 } else {
761 hlist[ipar] = new TH1D(name, title, nbins, bins->fArray);
762 }
763 hlist[ipar]->GetXaxis()->SetTitle(fYaxis.GetTitle());
764 if(arr != nullptr) {
765 (*arr)[ipar] = hlist[ipar];
766 }
767 }
768 snprintf(name, 2000, "%s_chi2", GetName());
769 delete gDirectory->FindObject(name);
770 TH1D* hchi2 = nullptr;
771 if(bins->fN == 0) {
772 hchi2 = new TH1D(name, "chisquare", nbins, fYaxis.GetXmin(), fYaxis.GetXmax());
773 } else {
774 hchi2 = new TH1D(name, "chisquare", nbins, bins->fArray);
775 }
776 hchi2->GetXaxis()->SetTitle(fYaxis.GetTitle());
777 if(arr != nullptr) {
778 (*arr)[npar] = hchi2;
779 }
780
781 // Loop on all bins in Y, generate a projection along X
782 // in case of sliding merge nstep=1, i.e. do slices starting for every bin
783 // now do not slices case with overflow (makes more sense)
784 for(Int_t bin = firstbin; bin + ngroup - 1 <= lastbin; bin += nstep) {
785 TH1D* proj = Projection("_temp", bin, bin + ngroup - 1, "e");
786 if(proj == nullptr) {
787 continue;
788 }
789 auto nentries = static_cast<Long64_t>(proj->GetEntries());
790 if(nentries == 0 || nentries < cut) {
791 delete proj;
792 continue;
793 }
794 f1->SetParameters(parsave);
795 proj->Fit(f1, opt.Data());
796 Int_t npfits = f1->GetNumberFitPoints();
797 if(npfits > npar && npfits >= cut) {
798 Int_t binOn = bin + ngroup / 2;
799 for(Int_t ipar = 0; ipar < npar; ++ipar) {
800 hlist[ipar]->Fill(fYaxis.GetBinCenter(binOn), f1->GetParameter(ipar));
801 hlist[ipar]->SetBinError(binOn, f1->GetParError(ipar));
802 }
803 hchi2->Fill(fYaxis.GetBinCenter(binOn), f1->GetChisquare() / (npfits - npar));
804 }
805 delete proj;
806 }
807 delete[] parsave;
808 delete[] name;
809 delete[] title;
810 delete[] hlist;
811}
812
813Int_t GHSym::GetBin(Int_t binx, Int_t biny, Int_t) const
814{
815 Int_t nBin = fXaxis.GetNbins() + 2;
816 if(binx < 0) {
817 binx = 0;
818 }
819 if(binx >= nBin) {
820 binx = nBin - 1;
821 }
822 if(biny < 0) {
823 biny = 0;
824 }
825 if(biny >= nBin) {
826 biny = nBin - 1;
827 }
828 if(biny <= binx) {
829 return binx + biny * (2 * fXaxis.GetNbins() - biny + 3) / 2;
830 }
831 return biny + binx * (2 * fXaxis.GetNbins() - binx + 3) / 2;
832}
833
834Double_t GHSym::GetBinWithContent2(Double_t content, Int_t& binx, Int_t& biny, Int_t firstxbin, Int_t lastxbin, Int_t firstybin, Int_t lastybin, Double_t maxdiff) const
835{
836 // compute first cell (binx,biny) in the range [firstxbin,lastxbin][firstybin,lastybin] for which
837 // diff = abs(cell_content-c) <= maxdiff
838 // In case several cells in the specified range with diff=0 are found
839 // the first cell found is returned in binx,biny.
840 // In case several cells in the specified range satisfy diff <=maxdiff
841 // the cell with the smallest difference is returned in binx,biny.
842 // In all cases the function returns the smallest difference.
843 //
844 // NOTE1: if firstxbin < 0, firstxbin is set to 1
845 // if(lastxbin < firstxbin then lastxbin is set to the number of bins in X
846 // ie if firstxbin=1 and lastxbin=0 (default) the search is on all bins in X except
847 // for X's under- and overflow bins.
848 // if firstybin < 0, firstybin is set to 1
849 // if(lastybin < firstybin then lastybin is set to the number of bins in Y
850 // ie if firstybin=1 and lastybin=0 (default) the search is on all bins in Y except
851 // for Y's under- and overflow bins.
852 // NOTE2: if maxdiff=0 (default), the first cell with content=c is returned.
853
854 if(fDimension != 2) {
855 binx = -1;
856 biny = -1;
857 Error("GetBinWithContent2", "function is only valid for 2-D histograms");
858 return 0;
859 }
860 if(firstxbin < 0) {
861 firstxbin = 1;
862 }
863 if(lastxbin < firstxbin) {
864 lastxbin = fXaxis.GetNbins();
865 }
866 if(firstybin < 0) {
867 firstybin = 1;
868 }
869 if(lastybin < firstybin) {
870 lastybin = fYaxis.GetNbins();
871 }
872 Double_t curmax = 1e240;
873 for(Int_t j = firstybin; j <= lastybin; j++) {
874 for(Int_t i = firstxbin; i <= lastxbin; i++) {
875 Double_t diff = TMath::Abs(GetBinContent(i, j) - content);
876 if(diff <= 0) {
877 binx = i;
878 biny = j;
879 return diff;
880 }
881 if(diff < curmax && diff <= maxdiff) {
882 curmax = diff, binx = i;
883 biny = j;
884 }
885 }
886 }
887 return curmax;
888}
889
890Double_t GHSym::GetCellContent(Int_t binx, Int_t biny) const
891{
892 return GetBinContent(GetBin(binx, biny));
893}
894
895Double_t GHSym::GetCellError(Int_t binx, Int_t biny) const
896{
897 return GetBinError(GetBin(binx, biny));
898}
899
900Double_t GHSym::GetCorrelationFactor(Int_t axis1, Int_t axis2) const
901{
902 //*-*-*-*-*-*-*-*Return correlation factor between axis1 and axis2*-*-*-*-*
903 //*-* ====================================================
904 if(axis1 < 1 || axis2 < 1 || axis1 > 2 || axis2 > 2) {
905 Error("GetCorrelationFactor", "Wrong parameters");
906 return 0;
907 }
908 if(axis1 == axis2) {
909 return 1;
910 }
911 Double_t rms1 = GetRMS(axis1);
912 if(rms1 == 0) {
913 return 0;
914 }
915 Double_t rms2 = GetRMS(axis2);
916 if(rms2 == 0) {
917 return 0;
918 }
919 return GetCovariance(axis1, axis2) / rms1 / rms2;
920}
921
922Double_t GHSym::GetCovariance(Int_t axis1, Int_t axis2) const
923{
924 //*-*-*-*-*-*-*-*Return covariance between axis1 and axis2*-*-*-*-*
925 //*-* ====================================================
926
927 if(axis1 < 1 || axis2 < 1 || axis1 > 2 || axis2 > 2) {
928 Error("GetCovariance", "Wrong parameters");
929 return 0;
930 }
931 std::array<Double_t, kNstat> stats;
932 GetStats(stats.data());
933 Double_t sumw = stats[0];
934 // Double_t sumw2 = stats[1];
935 Double_t sumwx = stats[2];
936 Double_t sumwx2 = stats[3];
937 Double_t sumwy = stats[4];
938 Double_t sumwy2 = stats[5];
939 Double_t sumwxy = stats[6];
940
941 if(sumw == 0) {
942 return 0;
943 }
944 if(axis1 == 1 && axis2 == 1) {
945 return TMath::Abs(sumwx2 / sumw - sumwx / sumw * sumwx / sumw);
946 }
947 if(axis1 == 2 && axis2 == 2) {
948 return TMath::Abs(sumwy2 / sumw - sumwy / sumw * sumwy / sumw);
949 }
950 return sumwxy / sumw - sumwx / sumw * sumwy / sumw;
951}
952
953void GHSym::GetRandom2(Double_t& x, Double_t& y)
954{
955 // return 2 random numbers along axis x and y distributed according
956 // the cellcontents of a 2-dim histogram
957 // return a NaN if the histogram has a bin with negative content
958
959 Int_t nbinsx = GetNbinsX();
960 Int_t nbinsy = GetNbinsY();
961 Int_t nbins = nbinsx * nbinsy;
962 Double_t integral = 0.;
963 // compute integral checking that all bins have positive content (see ROOT-5894)
964 if(fIntegral != nullptr) {
965 if(fIntegral[nbins + 1] != fEntries) {
966 integral = ComputeIntegral(true);
967 } else {
968 integral = fIntegral[nbins];
969 }
970 } else {
971 integral = ComputeIntegral(true);
972 }
973 if(integral == 0) {
974 x = 0;
975 y = 0;
976 return;
977 }
978 // case histogram has negative bins
979 if(integral == TMath::QuietNaN()) {
980 x = TMath::QuietNaN();
981 y = TMath::QuietNaN();
982 return;
983 }
984
985 Double_t r1 = gRandom->Rndm();
986 auto ibin = static_cast<int>(TMath::BinarySearch(nbins, fIntegral, r1));
987 Int_t biny = ibin / nbinsx;
988 Int_t binx = ibin - nbinsx * biny;
989 x = fXaxis.GetBinLowEdge(binx + 1);
990 if(r1 > fIntegral[ibin]) {
991 x += fXaxis.GetBinWidth(binx + 1) * (r1 - fIntegral[ibin]) / (fIntegral[ibin + 1] - fIntegral[ibin]);
992 }
993 y = fYaxis.GetBinLowEdge(biny + 1) + fYaxis.GetBinWidth(biny + 1) * gRandom->Rndm();
994}
995
996void GHSym::GetStats(Double_t* stats) const
997{
998 // fill the array stats from the contents of this histogram
999 // The array stats must be correctly dimensionned in the calling program.
1000 // stats[0] = sumw
1001 // stats[1] = sumw2
1002 // stats[2] = sumwx
1003 // stats[3] = sumwx2
1004 // stats[4] = sumwy
1005 // stats[5] = sumwy2
1006 // stats[6] = sumwxy
1007 //
1008 // If no axis-subranges are specified (via TAxis::SetRange), the array stats
1009 // is simply a copy of the statistics quantities computed at filling time.
1010 // If sub-ranges are specified, the function recomputes these quantities
1011 // from the bin contents in the current axis ranges.
1012 //
1013 // Note that the mean value/RMS is computed using the bins in the currently
1014 // defined ranges (see TAxis::SetRange). By default the ranges include
1015 // all bins from 1 to nbins included, excluding underflows and overflows.
1016 // To force the underflows and overflows in the computation, one must
1017 // call the static function TH1::StatOverflows(kTRUE) before filling
1018 // the histogram.
1019
1020 if(fBuffer != nullptr) {
1021 const_cast<GHSym*>(this)->BufferEmpty(); // NOLINT(cppcoreguidelines-pro-type-const-cast)
1022 }
1023
1024 if((fTsumw == 0 && fEntries > 0) || fXaxis.TestBit(TAxis::kAxisRange) || fYaxis.TestBit(TAxis::kAxisRange)) {
1025 for(Int_t bin = 0; bin < 7; ++bin) {
1026 stats[bin] = 0;
1027 }
1028
1029 Int_t firstBinX = fXaxis.GetFirst();
1030 Int_t lastBinX = fXaxis.GetLast();
1031 Int_t firstBinY = fYaxis.GetFirst();
1032 Int_t lastBinY = fYaxis.GetLast();
1033 // include underflow/overflow if TH1::StatOverflows(kTRUE) in case no range is set on the axis
1034 if(fgStatOverflows) {
1035 if(!fXaxis.TestBit(TAxis::kAxisRange)) {
1036 if(firstBinX == 1) {
1037 firstBinX = 0;
1038 }
1039 if(lastBinX == fXaxis.GetNbins()) {
1040 lastBinX += 1;
1041 }
1042 }
1043 if(!fYaxis.TestBit(TAxis::kAxisRange)) {
1044 if(firstBinY == 1) {
1045 firstBinY = 0;
1046 }
1047 if(lastBinY == fYaxis.GetNbins()) {
1048 lastBinY += 1;
1049 }
1050 }
1051 }
1052 for(Int_t biny = firstBinY; biny <= lastBinY; ++biny) {
1053 Double_t y = fYaxis.GetBinCenter(biny);
1054 for(Int_t binx = firstBinX; binx <= lastBinX; ++binx) {
1055 Int_t bin = GetBin(binx, biny);
1056 Double_t x = fXaxis.GetBinCenter(binx);
1057 Double_t w = GetBinContent(bin);
1058 Double_t err = TMath::Abs(GetBinError(bin));
1059 stats[0] += w;
1060 stats[1] += err * err;
1061 stats[2] += w * x;
1062 stats[3] += w * x * x;
1063 stats[4] += w * y;
1064 stats[5] += w * y * y;
1065 stats[6] += w * x * y;
1066 }
1067 }
1068 } else {
1069 stats[0] = fTsumw;
1070 stats[1] = fTsumw2;
1071 stats[2] = fTsumwx;
1072 stats[3] = fTsumwx2;
1073 stats[4] = fTsumwy;
1074 stats[5] = fTsumwy2;
1075 stats[6] = fTsumwxy;
1076 }
1077}
1078
1079Double_t GHSym::Integral(Option_t* option) const
1080{
1081 // Return integral of bin contents. Only bins in the bins range are considered.
1082 // By default the integral is computed as the sum of bin contents in the range.
1083 // if option "width" is specified, the integral is the sum of
1084 // the bin contents multiplied by the bin width in x and in y.
1085
1086 return Integral(fXaxis.GetFirst(), fXaxis.GetLast(), fYaxis.GetFirst(), fYaxis.GetLast(), option);
1087}
1088
1089Double_t GHSym::Integral(Int_t firstxbin, Int_t lastxbin, Int_t firstybin, Int_t lastybin, Option_t* option) const
1090{
1091 // Return integral of bin contents in range [firstxbin,lastxbin],[firstybin,lastybin]
1092 // for a 2-D histogram
1093 // By default the integral is computed as the sum of bin contents in the range.
1094 // if option "width" is specified, the integral is the sum of
1095 // the bin contents multiplied by the bin width in x and in y.
1096 double err = 0;
1097 return DoIntegral(firstxbin, lastxbin, firstybin, lastybin, err, option);
1098}
1099
1100Double_t GHSym::IntegralAndError(Int_t firstxbin, Int_t lastxbin, Int_t firstybin, Int_t lastybin, Double_t& error,
1101 Option_t* option) const
1102{
1103 // Return integral of bin contents in range [firstxbin,lastxbin],[firstybin,lastybin]
1104 // for a 2-D histogram. Calculates also the integral error using error propagation
1105 // from the bin errors assumming that all the bins are uncorrelated.
1106 // By default the integral is computed as the sum of bin contents in the range.
1107 // if option "width" is specified, the integral is the sum of
1108 // the bin contents multiplied by the bin width in x and in y.
1109
1110 return DoIntegral(firstxbin, lastxbin, firstybin, lastybin, error, option, kTRUE);
1111}
1112
1113#if ROOT_VERSION_CODE < ROOT_VERSION(6, 20, 0)
1114Double_t GHSym::Interpolate(Double_t)
1115#else
1116Double_t GHSym::Interpolate(Double_t) const
1117#endif
1118{
1119 // illegal for a TH2
1120 Error("Interpolate", "This function must be called with 2 arguments for a TH2");
1121 return 0;
1122}
1123
1124#if ROOT_VERSION_CODE < ROOT_VERSION(6, 20, 0)
1125Double_t GHSym::Interpolate(Double_t x, Double_t y)
1126#else
1127Double_t GHSym::Interpolate(Double_t x, Double_t y) const
1128#endif
1129{
1130 // Given a point P(x,y), Interpolate approximates the value via bilinear
1131 // interpolation based on the four nearest bin centers
1132 // see Wikipedia, Bilinear Interpolation
1133 // Andy Mastbaum 10/8/2008
1134 // vaguely based on R.Raja 6-Sep-2008
1135
1136 Double_t f = 0;
1137 Double_t x1 = 0;
1138 Double_t x2 = 0;
1139 Double_t y1 = 0;
1140 Double_t y2 = 0;
1141 Double_t dx = 0.;
1142 Double_t dy = 0.;
1143 Int_t bin_x = fXaxis.FindBin(x);
1144 Int_t bin_y = fYaxis.FindBin(y);
1145 if(bin_x < 1 || bin_x > GetNbinsX() || bin_y < 1 || bin_y > GetNbinsY()) {
1146 Error("Interpolate", "Cannot interpolate outside histogram domain.");
1147 return 0;
1148 }
1149 Int_t quadrant = 0; // CCW from UR 1,2,3,4
1150 // which quadrant of the bin (bin_P) are we in?
1151 dx = fXaxis.GetBinUpEdge(bin_x) - x;
1152 dy = fYaxis.GetBinUpEdge(bin_y) - y;
1153 if(dx <= fXaxis.GetBinWidth(bin_x) / 2 && dy <= fYaxis.GetBinWidth(bin_y) / 2) {
1154 quadrant = 1; // upper right
1155 }
1156 if(dx > fXaxis.GetBinWidth(bin_x) / 2 && dy <= fYaxis.GetBinWidth(bin_y) / 2) {
1157 quadrant = 2; // upper left
1158 }
1159 if(dx > fXaxis.GetBinWidth(bin_x) / 2 && dy > fYaxis.GetBinWidth(bin_y) / 2) {
1160 quadrant = 3; // lower left
1161 }
1162 if(dx <= fXaxis.GetBinWidth(bin_x) / 2 && dy > fYaxis.GetBinWidth(bin_y) / 2) {
1163 quadrant = 4; // lower right
1164 }
1165 switch(quadrant) {
1166 case 1:
1167 x1 = fXaxis.GetBinCenter(bin_x);
1168 y1 = fYaxis.GetBinCenter(bin_y);
1169 x2 = fXaxis.GetBinCenter(bin_x + 1);
1170 y2 = fYaxis.GetBinCenter(bin_y + 1);
1171 break;
1172 case 2:
1173 x1 = fXaxis.GetBinCenter(bin_x - 1);
1174 y1 = fYaxis.GetBinCenter(bin_y);
1175 x2 = fXaxis.GetBinCenter(bin_x);
1176 y2 = fYaxis.GetBinCenter(bin_y + 1);
1177 break;
1178 case 3:
1179 x1 = fXaxis.GetBinCenter(bin_x - 1);
1180 y1 = fYaxis.GetBinCenter(bin_y - 1);
1181 x2 = fXaxis.GetBinCenter(bin_x);
1182 y2 = fYaxis.GetBinCenter(bin_y);
1183 break;
1184 case 4:
1185 x1 = fXaxis.GetBinCenter(bin_x);
1186 y1 = fYaxis.GetBinCenter(bin_y - 1);
1187 x2 = fXaxis.GetBinCenter(bin_x + 1);
1188 y2 = fYaxis.GetBinCenter(bin_y);
1189 break;
1190 }
1191 Int_t bin_x1 = fXaxis.FindBin(x1);
1192 if(bin_x1 < 1) {
1193 bin_x1 = 1;
1194 }
1195 Int_t bin_x2 = fXaxis.FindBin(x2);
1196 if(bin_x2 > GetNbinsX()) {
1197 bin_x2 = GetNbinsX();
1198 }
1199 Int_t bin_y1 = fYaxis.FindBin(y1);
1200 if(bin_y1 < 1) {
1201 bin_y1 = 1;
1202 }
1203 Int_t bin_y2 = fYaxis.FindBin(y2);
1204 if(bin_y2 > GetNbinsY()) {
1205 bin_y2 = GetNbinsY();
1206 }
1207 Int_t bin_q22 = GetBin(bin_x2, bin_y2);
1208 Int_t bin_q12 = GetBin(bin_x1, bin_y2);
1209 Int_t bin_q11 = GetBin(bin_x1, bin_y1);
1210 Int_t bin_q21 = GetBin(bin_x2, bin_y1);
1211 Double_t q11 = GetBinContent(bin_q11);
1212 Double_t q12 = GetBinContent(bin_q12);
1213 Double_t q21 = GetBinContent(bin_q21);
1214 Double_t q22 = GetBinContent(bin_q22);
1215 Double_t d = 1.0 * (x2 - x1) * (y2 - y1);
1216 f = 1.0 * q11 / d * (x2 - x) * (y2 - y) + 1.0 * q21 / d * (x - x1) * (y2 - y) + 1.0 * q12 / d * (x2 - x) * (y - y1) +
1217 1.0 * q22 / d * (x - x1) * (y - y1);
1218 return f;
1219}
1220
1221#if ROOT_VERSION_CODE < ROOT_VERSION(6, 20, 0)
1222Double_t GHSym::Interpolate(Double_t, Double_t, Double_t)
1223#else
1224Double_t GHSym::Interpolate(Double_t, Double_t, Double_t) const
1225#endif
1226{
1227 // illegal for a TH2
1228 Error("Interpolate", "This function must be called with 2 arguments for a TH2");
1229 return 0;
1230}
1231
1232Double_t GHSym::KolmogorovTest(const TH1* h2, Option_t* option) const
1233{
1234 // Statistical test of compatibility in shape between
1235 // THIS histogram and h2, using Kolmogorov test.
1236 // Default: Ignore under- and overflow bins in comparison
1237 //
1238 // option is a character string to specify options
1239 // "U" include Underflows in test
1240 // "O" include Overflows
1241 // "N" include comparison of normalizations
1242 // "D" Put out a line of "Debug" printout
1243 // "M" Return the Maximum Kolmogorov distance instead of prob
1244 //
1245 // The returned function value is the probability of test
1246 // (much less than one means NOT compatible)
1247 //
1248 // The KS test uses the distance between the pseudo-CDF's obtained
1249 // from the histogram. Since in 2D the order for generating the pseudo-CDF is
1250 // arbitrary, two pairs of pseudo-CDF are used, one starting from the x axis the
1251 // other from the y axis and the maximum distance is the average of the two maximum
1252 // distances obtained.
1253 //
1254 // Code adapted by Rene Brun from original HBOOK routine HDIFF
1255
1256 TString opt = option;
1257 opt.ToUpper();
1258
1259 Double_t prb = 0;
1260 TH1* h1 = const_cast<TH1*>(static_cast<const TH1*>(this)); // NOLINT(cppcoreguidelines-pro-type-const-cast)
1261 if(h2 == nullptr) {
1262 return 0;
1263 }
1264 TAxis* xaxis1 = h1->GetXaxis();
1265 auto* xaxis2 = const_cast<TAxis*>(h2->GetXaxis()); // NOLINT(cppcoreguidelines-pro-type-const-cast)
1266 TAxis* yaxis1 = h1->GetYaxis();
1267 auto* yaxis2 = const_cast<TAxis*>(h2->GetYaxis()); // NOLINT(cppcoreguidelines-pro-type-const-cast)
1268 Int_t ncx1 = xaxis1->GetNbins();
1269 Int_t ncx2 = xaxis2->GetNbins();
1270 Int_t ncy1 = yaxis1->GetNbins();
1271 Int_t ncy2 = yaxis2->GetNbins();
1272
1273 // Check consistency of dimensions
1274 if(h1->GetDimension() != 2 || h2->GetDimension() != 2) {
1275 Error("KolmogorovTest", "Histograms must be 2-D\n");
1276 return 0;
1277 }
1278
1279 // Check consistency in number of channels
1280 if(ncx1 != ncx2) {
1281 Error("KolmogorovTest", "Number of channels in X is different, %d and %d\n", ncx1, ncx2);
1282 return 0;
1283 }
1284 if(ncy1 != ncy2) {
1285 Error("KolmogorovTest", "Number of channels in Y is different, %d and %d\n", ncy1, ncy2);
1286 return 0;
1287 }
1288
1289 // Check consistency in channel edges
1290 Bool_t afunc1 = kFALSE;
1291 Bool_t afunc2 = kFALSE;
1292 Double_t difprec = 1e-5;
1293 Double_t diff1 = TMath::Abs(xaxis1->GetXmin() - xaxis2->GetXmin());
1294 Double_t diff2 = TMath::Abs(xaxis1->GetXmax() - xaxis2->GetXmax());
1295 if(diff1 > difprec || diff2 > difprec) {
1296 Error("KolmogorovTest", "histograms with different binning along X");
1297 return 0;
1298 }
1299 diff1 = TMath::Abs(yaxis1->GetXmin() - yaxis2->GetXmin());
1300 diff2 = TMath::Abs(yaxis1->GetXmax() - yaxis2->GetXmax());
1301 if(diff1 > difprec || diff2 > difprec) {
1302 Error("KolmogorovTest", "histograms with different binning along Y");
1303 return 0;
1304 }
1305
1306 // Should we include Uflows, Oflows?
1307 Int_t ibeg = 1;
1308 Int_t jbeg = 1;
1309 Int_t iend = ncx1;
1310 Int_t jend = ncy1;
1311 if(opt.Contains("U")) {
1312 ibeg = 0;
1313 jbeg = 0;
1314 }
1315 if(opt.Contains("O")) {
1316 iend = ncx1 + 1;
1317 jend = ncy1 + 1;
1318 }
1319
1320 Double_t sum1 = 0;
1321 Double_t sum2 = 0;
1322 Double_t w1 = 0;
1323 Double_t w2 = 0;
1324 for(Int_t i = ibeg; i <= iend; ++i) {
1325 for(Int_t j = jbeg; j <= jend; ++j) {
1326 sum1 += h1->GetBinContent(i, j);
1327 sum2 += h2->GetBinContent(i, j);
1328 Double_t ew1 = h1->GetBinError(i, j);
1329 Double_t ew2 = h2->GetBinError(i, j);
1330 w1 += ew1 * ew1;
1331 w2 += ew2 * ew2;
1332 }
1333 }
1334
1335 // Check that both scatterplots contain events
1336 if(sum1 == 0) {
1337 Error("KolmogorovTest", "Integral is zero for h1=%s\n", h1->GetName());
1338 return 0;
1339 }
1340 if(sum2 == 0) {
1341 Error("KolmogorovTest", "Integral is zero for h2=%s\n", h2->GetName());
1342 return 0;
1343 }
1344 // calculate the effective entries.
1345 // the case when errors are zero (w1 == 0 or w2 ==0) are equivalent to
1346 // compare to a function. In that case the rescaling is done only on sqrt(esum2) or sqrt(esum1)
1347 Double_t esum1 = 0.;
1348 Double_t esum2 = 0.;
1349 if(w1 > 0) {
1350 esum1 = sum1 * sum1 / w1;
1351 } else {
1352 afunc1 = kTRUE; // use later for calculating z
1353 }
1354
1355 if(w2 > 0) {
1356 esum2 = sum2 * sum2 / w2;
1357 } else {
1358 afunc2 = kTRUE; // use later for calculating z
1359 }
1360
1361 if(afunc2 && afunc1) {
1362 Error("KolmogorovTest", "Errors are zero for both histograms\n");
1363 return 0;
1364 }
1365
1366 // Find first Kolmogorov distance
1367 Double_t s1 = 1 / sum1;
1368 Double_t s2 = 1 / sum2;
1369 Double_t dfmax1 = 0;
1370 Double_t rsum1 = 0.;
1371 Double_t rsum2 = 0.;
1372 for(Int_t i = ibeg; i <= iend; ++i) {
1373 for(Int_t j = jbeg; j <= jend; ++j) {
1374 rsum1 += s1 * h1->GetCellContent(i, j);
1375 rsum2 += s2 * h2->GetCellContent(i, j);
1376 dfmax1 = TMath::Max(dfmax1, TMath::Abs(rsum1 - rsum2));
1377 }
1378 }
1379
1380 // Find second Kolmogorov distance
1381 Double_t dfmax2 = 0;
1382 rsum1 = 0.;
1383 rsum2 = 0.;
1384 for(Int_t j = jbeg; j <= jend; ++j) {
1385 for(Int_t i = ibeg; i <= iend; ++i) {
1386 rsum1 += s1 * h1->GetCellContent(i, j);
1387 rsum2 += s2 * h2->GetCellContent(i, j);
1388 dfmax2 = TMath::Max(dfmax2, TMath::Abs(rsum1 - rsum2));
1389 }
1390 }
1391
1392 // Get Kolmogorov probability: use effective entries, esum1 or esum2, for normalizing it
1393 Double_t factnm = 0.;
1394 if(afunc1) {
1395 factnm = TMath::Sqrt(esum2);
1396 } else if(afunc2) {
1397 factnm = TMath::Sqrt(esum1);
1398 } else {
1399 factnm = TMath::Sqrt(esum1 * sum2 / (esum1 + esum2));
1400 }
1401
1402 // take average of the two distances
1403 Double_t dfmax = 0.5 * (dfmax1 + dfmax2);
1404 Double_t z = dfmax * factnm;
1405
1406 prb = TMath::KolmogorovProb(z);
1407
1408 Double_t prb1 = 0;
1409 Double_t prb2 = 0;
1410 // option N to combine normalization makes sense if both afunc1 and afunc2 are false
1411 if(opt.Contains("N") && !(afunc1 || afunc2)) {
1412 // Combine probabilities for shape and normalization
1413 prb1 = prb;
1414 Double_t d12 = esum1 - esum2;
1415 Double_t chi2 = d12 * d12 / (esum1 + esum2);
1416 prb2 = TMath::Prob(chi2, 1);
1417 // see Eadie et al., section 11.6.2
1418 if(prb > 0 && prb2 > 0) {
1419 prb = prb * prb2 * (1 - TMath::Log(prb * prb2));
1420 } else {
1421 prb = 0;
1422 }
1423 }
1424 // debug printout
1425 if(opt.Contains("D")) {
1426 std::cout << " Kolmo Prob h1 = " << h1->GetName() << ", sum1 = " << sum1 << std::endl;
1427 std::cout << " Kolmo Prob h2 = " << h2->GetName() << ", sum2 = " << sum2 << std::endl;
1428 std::cout << " Kolmo Probabil = " << prb << ", Max dist = " << dfmax << std::endl;
1429 if(opt.Contains("N")) {
1430 std::cout << " Kolmo Probabil = " << prb1 << " for shape alone, " << prb2 << " for normalisation alone" << std::endl;
1431 }
1432 }
1433 // This numerical error condition should never occur:
1434 if(TMath::Abs(rsum1 - 1) > 0.002) {
1435 Warning("KolmogorovTest", "Numerical problems with h1=%s\n", h1->GetName());
1436 }
1437 if(TMath::Abs(rsum2 - 1) > 0.002) {
1438 Warning("KolmogorovTest", "Numerical problems with h2=%s\n", h2->GetName());
1439 }
1440
1441 if(opt.Contains("M")) {
1442 return dfmax; // return avergae of max distance
1443 }
1444
1445 return prb;
1446}
1447
1448Long64_t GHSym::Merge(TCollection* list)
1449{
1450 // Add all histograms in the collection to this histogram.
1451 // This function computes the min/max for the axes,
1452 // compute a new number of bins, if necessary,
1453 // add bin contents, errors and statistics.
1454 // If overflows are present and limits are different the function will fail.
1455 // The function returns the total number of entries in the result histogram
1456 // if the merge is successfull, -1 otherwise.
1457 //
1458 // IMPORTANT remark. The 2 axis x and y may have different number
1459 // of bins and different limits, BUT the largest bin width must be
1460 // a multiple of the smallest bin width and the upper limit must also
1461 // be a multiple of the bin width.
1462
1463 if(list == nullptr) {
1464 return 0;
1465 }
1466 if(list->IsEmpty()) {
1467 return static_cast<Long64_t>(GetEntries());
1468 }
1469
1470 TList inlist;
1471 inlist.AddAll(list);
1472
1473 TAxis newXAxis;
1474 TAxis newYAxis;
1475 Bool_t initialLimitsFound = kFALSE;
1476 Bool_t allSameLimits = kTRUE;
1477 Bool_t sameLimitsX = kTRUE;
1478 Bool_t sameLimitsY = kTRUE;
1479 Bool_t allHaveLimits = kTRUE;
1480 Bool_t firstHistWithLimits = kTRUE;
1481
1482 TIter next(&inlist);
1483 GHSym* h = this;
1484 do {
1485 Bool_t hasLimits = h->GetXaxis()->GetXmin() < h->GetXaxis()->GetXmax();
1486 allHaveLimits = allHaveLimits && hasLimits;
1487
1488 if(hasLimits) {
1489 h->BufferEmpty();
1490
1491 // this is done in case the first histograms are empty and
1492 // the histogram have different limits
1493 if(firstHistWithLimits) {
1494 // set axis limits in the case the first histogram did not have limits
1495 if(h != this) {
1496 if(!SameLimitsAndNBins(fXaxis, *(h->GetXaxis()))) {
1497 if(h->GetXaxis()->GetXbins()->GetSize() != 0) {
1498 fXaxis.Set(h->GetXaxis()->GetNbins(), h->GetXaxis()->GetXbins()->GetArray());
1499 } else {
1500 fXaxis.Set(h->GetXaxis()->GetNbins(), h->GetXaxis()->GetXmin(), h->GetXaxis()->GetXmax());
1501 }
1502 }
1503 if(!SameLimitsAndNBins(fYaxis, *(h->GetYaxis()))) {
1504 if(h->GetYaxis()->GetXbins()->GetSize() != 0) {
1505 fYaxis.Set(h->GetYaxis()->GetNbins(), h->GetYaxis()->GetXbins()->GetArray());
1506 } else {
1507 fYaxis.Set(h->GetYaxis()->GetNbins(), h->GetYaxis()->GetXmin(), h->GetYaxis()->GetXmax());
1508 }
1509 }
1510 }
1511 firstHistWithLimits = kFALSE;
1512 }
1513
1514 if(!initialLimitsFound) {
1515 // this is executed the first time an histogram with limits is found
1516 // to set some initial values on the new axes
1517 initialLimitsFound = kTRUE;
1518 if(h->GetXaxis()->GetXbins()->GetSize() != 0) {
1519 newXAxis.Set(h->GetXaxis()->GetNbins(), h->GetXaxis()->GetXbins()->GetArray());
1520 } else {
1521 newXAxis.Set(h->GetXaxis()->GetNbins(), h->GetXaxis()->GetXmin(), h->GetXaxis()->GetXmax());
1522 }
1523 if(h->GetYaxis()->GetXbins()->GetSize() != 0) {
1524 newYAxis.Set(h->GetYaxis()->GetNbins(), h->GetYaxis()->GetXbins()->GetArray());
1525 } else {
1526 newYAxis.Set(h->GetYaxis()->GetNbins(), h->GetYaxis()->GetXmin(), h->GetYaxis()->GetXmax());
1527 }
1528 } else {
1529 // check first if histograms have same bins in X
1530 if(!SameLimitsAndNBins(newXAxis, *(h->GetXaxis()))) {
1531 sameLimitsX = kFALSE;
1532 // recompute in this case the optimal limits
1533 // The condition to works is that the histogram have same bin with
1534 // and one common bin edge
1535 if(!RecomputeAxisLimits(newXAxis, *(h->GetXaxis()))) {
1536 Error("Merge", "Cannot merge histograms - limits are inconsistent:\n "
1537 "first: (%d, %f, %f), second: (%d, %f, %f)",
1538 newXAxis.GetNbins(), newXAxis.GetXmin(), newXAxis.GetXmax(), h->GetXaxis()->GetNbins(),
1539 h->GetXaxis()->GetXmin(), h->GetXaxis()->GetXmax());
1540 return -1;
1541 }
1542 }
1543
1544 // check first if histograms have same bins in Y
1545 if(!SameLimitsAndNBins(newYAxis, *(h->GetYaxis()))) {
1546 sameLimitsY = kFALSE;
1547 // recompute in this case the optimal limits
1548 // The condition to works is that the histogram have same bin with
1549 // and one common bin edge
1550 if(!RecomputeAxisLimits(newYAxis, *(h->GetYaxis()))) {
1551 Error("Merge", "Cannot merge histograms - limits are inconsistent:\n "
1552 "first: (%d, %f, %f), second: (%d, %f, %f)",
1553 newYAxis.GetNbins(), newYAxis.GetXmin(), newYAxis.GetXmax(), h->GetYaxis()->GetNbins(),
1554 h->GetYaxis()->GetXmin(), h->GetYaxis()->GetXmax());
1555 return -1;
1556 }
1557 }
1558 allSameLimits = sameLimitsY && sameLimitsX;
1559 }
1560 }
1561 } while((h = static_cast<GHSym*>(next())) != nullptr);
1562 if(h == nullptr && ((*next) != nullptr)) {
1563 Error("Merge", "Attempt to merge object of class: %s to a %s", (*next)->ClassName(), ClassName());
1564 return -1;
1565 }
1566 next.Reset();
1567
1568 // In the case of histogram with different limits
1569 // newX(Y)Axis will now have the new found limits
1570 // but one needs first to clone this histogram to perform the merge
1571 // The clone is not needed when all histograms have the same limits
1572 TH2* hclone = nullptr;
1573 if(!allSameLimits) {
1574 // We don't want to add the clone to gDirectory,
1575 // so remove our kMustCleanup bit temporarily
1576 Bool_t mustCleanup = TestBit(kMustCleanup);
1577 if(mustCleanup) {
1578 ResetBit(kMustCleanup);
1579 }
1580 hclone = static_cast<TH2*>(IsA()->New());
1581 hclone->SetDirectory(nullptr);
1582 Copy(*hclone);
1583 if(mustCleanup) {
1584 SetBit(kMustCleanup);
1585 }
1586 BufferEmpty(1); // To remove buffer.
1587 Reset(); // BufferEmpty sets limits so we can't use it later.
1588 SetEntries(0);
1589 inlist.AddFirst(hclone);
1590 }
1591
1592 if(!allSameLimits && initialLimitsFound) {
1593 if(!sameLimitsX) {
1594 fXaxis.SetRange(0, 0);
1595 if(newXAxis.GetXbins()->GetSize() != 0) {
1596 fXaxis.Set(newXAxis.GetNbins(), newXAxis.GetXbins()->GetArray());
1597 } else {
1598 fXaxis.Set(newXAxis.GetNbins(), newXAxis.GetXmin(), newXAxis.GetXmax());
1599 }
1600 }
1601 if(!sameLimitsY) {
1602 fYaxis.SetRange(0, 0);
1603 if(newYAxis.GetXbins()->GetSize() != 0) {
1604 fYaxis.Set(newYAxis.GetNbins(), newYAxis.GetXbins()->GetArray());
1605 } else {
1606 fYaxis.Set(newYAxis.GetNbins(), newYAxis.GetXmin(), newYAxis.GetXmax());
1607 }
1608 }
1609 fZaxis.Set(1, 0, 1);
1610 fNcells = (fXaxis.GetNbins() + 2) * (fYaxis.GetNbins() + 2);
1611 SetBinsLength(fNcells);
1612 if(fSumw2.fN != 0) {
1613 fSumw2.Set(fNcells);
1614 }
1615 }
1616
1617 if(!allHaveLimits) {
1618 // fill this histogram with all the data from buffers of histograms without limits
1619 while((h = static_cast<GHSym*>(next())) != nullptr) {
1620 if(h->GetXaxis()->GetXmin() >= h->GetXaxis()->GetXmax() && (h->fBuffer != nullptr)) {
1621 // no limits
1622 auto nbentries = static_cast<Int_t>(h->fBuffer[0]);
1623 for(Int_t i = 0; i < nbentries; i++) {
1624 Fill(h->fBuffer[3 * i + 2], h->fBuffer[3 * i + 3], h->fBuffer[3 * i + 1]);
1625 }
1626 // Entries from buffers have to be filled one by one
1627 // because FillN doesn't resize histograms.
1628 }
1629 }
1630 if(!initialLimitsFound) {
1631 if(hclone != nullptr) {
1632 inlist.Remove(hclone);
1633 delete hclone;
1634 }
1635 return static_cast<Long64_t>(GetEntries()); // all histograms have been processed
1636 }
1637 next.Reset();
1638 }
1639
1640 // merge bin contents and errors
1641 std::array<Double_t, kNstat> stats;
1642 std::array<Double_t, kNstat> totstats;
1643 for(Int_t i = 0; i < kNstat; ++i) {
1644 totstats[i] = stats[i] = 0;
1645 }
1646 GetStats(totstats.data());
1647 Double_t nentries = GetEntries();
1648 Int_t ix = 0;
1649 Int_t iy = 0;
1650 Double_t cu = 0.;
1651 Bool_t canExtend = CanExtendAllAxes();
1652 SetCanExtend(TH1::kNoAxis); // reset, otherwise setting the under/overflow will extend the axis
1653
1654 while((h = static_cast<GHSym*>(next())) != nullptr) {
1655 // skip empty histograms
1656 Double_t histEntries = h->GetEntries();
1657 if(h->fTsumw == 0 && histEntries == 0) {
1658 continue;
1659 }
1660
1661 // process only if the histogram has limits; otherwise it was processed before
1662 if(h->GetXaxis()->GetXmin() < h->GetXaxis()->GetXmax()) {
1663 // import statistics
1664 h->GetStats(stats.data());
1665 for(Int_t i = 0; i < kNstat; ++i) {
1666 totstats[i] += stats[i];
1667 }
1668 nentries += histEntries;
1669
1670 Int_t nx = h->GetXaxis()->GetNbins();
1671 Int_t ny = h->GetYaxis()->GetNbins();
1672
1673 for(Int_t biny = 0; biny <= ny + 1; ++biny) {
1674 if(!allSameLimits) {
1675 iy = fYaxis.FindBin(h->GetYaxis()->GetBinCenter(biny));
1676 } else {
1677 iy = biny;
1678 }
1679 for(Int_t binx = 0; binx <= nx + 1; ++binx) {
1680 cu = h->GetBinContent(binx, biny);
1681 if(!allSameLimits) {
1682 if(cu != 0 && ((!sameLimitsX && (binx == 0 || binx == nx + 1)) ||
1683 (!sameLimitsY && (biny == 0 || biny == ny + 1)))) {
1684 Error("Merge", "Cannot merge histograms - the histograms have"
1685 " different limits and undeflows/overflows are present."
1686 " The initial histogram is now broken!");
1687 return -1;
1688 }
1689 ix = fXaxis.FindBin(h->GetXaxis()->GetBinCenter(binx));
1690 } else {
1691 // case histograms with the same limits
1692 ix = binx;
1693 }
1694 Int_t ibin = GetBin(ix, iy);
1695
1696 if(ibin < 0) {
1697 continue;
1698 }
1699 AddBinContent(ibin, cu);
1700 if(fSumw2.fN != 0) {
1701 Double_t error1 = h->GetBinError(GetBin(binx, biny));
1702 fSumw2.fArray[ibin] += error1 * error1;
1703 }
1704 }
1705 }
1706 }
1707 }
1708 SetCanExtend(static_cast<UInt_t>(canExtend));
1709
1710 // copy merged stats
1711 PutStats(totstats.data());
1712 SetEntries(nentries);
1713 if(hclone != nullptr) {
1714 inlist.Remove(hclone);
1715 delete hclone;
1716 }
1717 return static_cast<Long64_t>(nentries);
1718}
1719
1720TProfile* GHSym::Profile(const char* name, Int_t firstbin, Int_t lastbin, Option_t* option) const
1721{
1722 // *-*-*-*-*Project a 2-D histogram into a profile histogram along X*-*-*-*-*-*
1723 // *-* ========================================================
1724 //
1725 // The projection is made from the channels along the Y axis
1726 // ranging from firstybin to lastybin included.
1727 // By default, bins 1 to ny are included
1728 // When all bins are included, the number of entries in the projection
1729 // is set to the number of entries of the 2-D histogram, otherwise
1730 // the number of entries is incremented by 1 for all non empty cells.
1731 //
1732 // if option "d" is specified, the profile is drawn in the current pad.
1733 //
1734 // if option "o" original axis range of the target axes will be
1735 // kept, but only bins inside the selected range will be filled.
1736 //
1737 // The option can also be used to specify the projected profile error type.
1738 // Values which can be used are 's', 'i', or 'g'. See TProfile::BuildOptions for details
1739 //
1740 // Using a TCutG object, it is possible to select a sub-range of a 2-D histogram.
1741 // One must create a graphical cut (mouse or C++) and specify the name
1742 // of the cut between [] in the option.
1743 // For example, with a TCutG named "cutg", one can call:
1744 // myhist->ProfileX(" ",firstybin,lastybin,"[cutg]");
1745 // To invert the cut, it is enough to put a "-" in front of its name:
1746 // myhist->ProfileX(" ",firstybin,lastybin,"[-cutg]");
1747 // It is possible to apply several cuts ("," means logical AND):
1748 // myhist->ProfileX(" ",firstybin,lastybin,"[cutg1,cutg2]");
1749 //
1750 // NOTE that if a TProfile named "name" exists in the current directory or pad with
1751 // a compatible axis the profile is reset and filled again with the projected contents of the TH2.
1752 // In the case of axis incompatibility an error is reported and a NULL pointer is returned.
1753 //
1754 // NOTE that the X axis attributes of the TH2 are copied to the X axis of the profile.
1755 //
1756 // NOTE that the default under- / overflow behavior differs from what ProjectionX
1757 // does! Profiles take the bin center into account, so here the under- and overflow
1758 // bins are ignored by default.
1759
1760 TString opt = option;
1761 // extract cut infor
1762 TString cut;
1763 Int_t i1 = opt.Index("[");
1764 if(i1 >= 0) {
1765 Int_t i2 = opt.Index("]");
1766 cut = opt(i1, i2 - i1 + 1);
1767 }
1768 opt.ToLower();
1769 bool originalRange = opt.Contains("o");
1770
1771 Int_t inN = fYaxis.GetNbins();
1772 const char* expectedName = "_pf";
1773
1774 Int_t firstOutBin = fXaxis.GetFirst();
1775 Int_t lastOutBin = fXaxis.GetLast();
1776 if(firstOutBin == 0 && lastOutBin == 0) {
1777 firstOutBin = 1;
1778 lastOutBin = fXaxis.GetNbins();
1779 }
1780
1781 if(lastbin < firstbin && fYaxis.TestBit(TAxis::kAxisRange)) {
1782 firstbin = fYaxis.GetFirst();
1783 lastbin = fYaxis.GetLast();
1784 // For special case of TAxis::SetRange, when first == 1 and last
1785 // = N and the range bit has been set, the TAxis will return 0
1786 // for both.
1787 if(firstbin == 0 && lastbin == 0) {
1788 firstbin = 1;
1789 lastbin = fYaxis.GetNbins();
1790 }
1791 }
1792 if(firstbin < 0) {
1793 firstbin = 1;
1794 }
1795 if(lastbin < 0) {
1796 lastbin = inN;
1797 }
1798 if(lastbin > inN + 1) {
1799 lastbin = inN;
1800 }
1801
1802 // Create the profile histogram
1803 char* pname = const_cast<char*>(name); // NOLINT(cppcoreguidelines-pro-type-const-cast)
1804 if((name != nullptr) && strcmp(name, expectedName) == 0) {
1805 auto nch = strlen(GetName()) + 5;
1806 pname = new char[nch];
1807 snprintf(pname, nch, "%s%s", GetName(), name);
1808 }
1809 TProfile* h1 = nullptr;
1810 // check if a profile with identical name exist
1811 // if compatible reset and re-use previous histogram
1812 TObject* h1obj = gROOT->FindObject(pname);
1813 if((h1obj != nullptr) && h1obj->InheritsFrom(TH1::Class())) {
1814 if(h1obj->IsA() != TProfile::Class()) {
1815 Error("DoProfile", "Histogram with name %s must be a TProfile and is a %s", name, h1obj->ClassName());
1816 return nullptr;
1817 }
1818 h1 = static_cast<TProfile*>(h1obj);
1819 // reset the existing histogram and set always the new binning for the axis
1820 // This avoid problems when the histogram already exists and the histograms is rebinned or its range has changed
1821 // (see https://savannah.cern.ch/bugs/?94101 or https://savannah.cern.ch/bugs/?95808 )
1822 h1->Reset();
1823 const TArrayD* xbins = fXaxis.GetXbins();
1824 if(xbins->fN == 0) {
1825 if(originalRange) {
1826 h1->SetBins(fXaxis.GetNbins(), fXaxis.GetXmin(), fXaxis.GetXmax());
1827 } else {
1828 h1->SetBins(lastOutBin - firstOutBin + 1, fXaxis.GetBinLowEdge(firstOutBin),
1829 fXaxis.GetBinUpEdge(lastOutBin));
1830 }
1831 } else {
1832 // case variable bins
1833 if(originalRange) {
1834 h1->SetBins(fXaxis.GetNbins(), xbins->fArray);
1835 } else {
1836 h1->SetBins(lastOutBin - firstOutBin + 1, &xbins->fArray[firstOutBin - 1]);
1837 }
1838 }
1839 }
1840
1841 Int_t ncuts = 0;
1842 if(opt.Contains("[")) {
1843 const_cast<GHSym*>(this)->GetPainter(); // NOLINT(cppcoreguidelines-pro-type-const-cast)
1844 if(fPainter != nullptr) {
1845 ncuts = fPainter->MakeCuts(const_cast<char*>(cut.Data())); // NOLINT(cppcoreguidelines-pro-type-const-cast)
1846 }
1847 }
1848
1849 if(h1 == nullptr) {
1850 const TArrayD* bins = fXaxis.GetXbins();
1851 if(bins->fN == 0) {
1852 if(originalRange) {
1853 h1 = new TProfile(pname, GetTitle(), fXaxis.GetNbins(), fXaxis.GetXmin(), fXaxis.GetXmax(), opt);
1854 } else {
1855 h1 = new TProfile(pname, GetTitle(), lastOutBin - firstOutBin + 1, fXaxis.GetBinLowEdge(firstOutBin),
1856 fXaxis.GetBinUpEdge(lastOutBin), opt);
1857 }
1858 } else {
1859 // case variable bins
1860 if(originalRange) {
1861 h1 = new TProfile(pname, GetTitle(), fXaxis.GetNbins(), bins->fArray, opt);
1862 } else {
1863 h1 = new TProfile(pname, GetTitle(), lastOutBin - firstOutBin + 1, &bins->fArray[firstOutBin - 1], opt);
1864 }
1865 }
1866 }
1867 if(pname != name) {
1868 delete[] pname;
1869 }
1870
1871 // Copy attributes
1872 h1->GetXaxis()->ImportAttributes(&fXaxis);
1873 h1->SetLineColor(GetLineColor());
1874 h1->SetFillColor(GetFillColor());
1875 h1->SetMarkerColor(GetMarkerColor());
1876 h1->SetMarkerStyle(GetMarkerStyle());
1877
1878 // check if histogram is weighted
1879 // in case need to store sum of weight square/bin for the profile
1880 bool useWeights = (GetSumw2N() > 0);
1881 if(useWeights) {
1882 h1->Sumw2();
1883 }
1884
1885 // Fill the profile histogram
1886 // no entries/bin is available so can fill only using bin content as weight
1887 TArrayD& binSumw2 = *(h1->GetBinSumw2());
1888
1889 // implement filling of projected histogram
1890 // outbin is bin number of fXaxis (the projected axis). Loop is done on all bin of TH2 histograms
1891 // inbin is the axis being integrated. Loop is done only on the selected bins
1892 for(Int_t outbin = 0; outbin <= fXaxis.GetNbins() + 1; ++outbin) {
1893 if(fXaxis.TestBit(TAxis::kAxisRange) && (outbin < firstOutBin || outbin > lastOutBin)) {
1894 continue;
1895 }
1896
1897 // find corresponding bin number in h1 for outbin (binOut)
1898 Double_t xOut = fXaxis.GetBinCenter(outbin);
1899 Int_t binOut = h1->GetXaxis()->FindBin(xOut);
1900 if(binOut < 0) {
1901 continue;
1902 }
1903
1904 for(Int_t inbin = firstbin; inbin <= lastbin; ++inbin) {
1905 Int_t binx = outbin;
1906 Int_t biny = inbin;
1907
1908 if(ncuts != 0) {
1909 if(!fPainter->IsInside(binx, biny)) {
1910 continue;
1911 }
1912 }
1913 Int_t bin = GetBin(binx, biny);
1914 Double_t cxy = GetBinContent(bin);
1915
1916 if(cxy != 0.0) {
1917 Double_t tmp = 0;
1918 // the following fill update wrongly the fBinSumw2- need to save it before
1919 if(useWeights) {
1920 tmp = binSumw2.fArray[binOut];
1921 }
1922 h1->Fill(xOut, fYaxis.GetBinCenter(inbin), cxy);
1923 if(useWeights) {
1924 binSumw2.fArray[binOut] = tmp + fSumw2.fArray[bin];
1925 }
1926 }
1927 }
1928 }
1929
1930 // the statistics must be recalculated since by using the Fill method the total sum of weight^2 is
1931 // not computed correctly
1932 // for a profile does not much sense to re-use statistics of original TH2
1933 h1->ResetStats();
1934 // Also we need to set the entries since they have not been correctly calculated during the projection
1935 // we can only set them to the effective entries
1936 h1->SetEntries(h1->GetEffectiveEntries());
1937
1938 if(opt.Contains("d")) {
1939 TVirtualPad* padsav = gPad;
1940 TVirtualPad* pad = gROOT->GetSelectedPad();
1941 if(pad != nullptr) {
1942 pad->cd();
1943 }
1944 opt.Remove(opt.First("d"), 1);
1945 if(!gPad || !gPad->FindObject(h1)) {
1946 h1->Draw(opt);
1947 } else {
1948 h1->Paint(opt);
1949 }
1950 if(padsav != nullptr) {
1951 padsav->cd();
1952 }
1953 }
1954 return h1;
1955}
1956
1957TH1D* GHSym::Projection(const char* name, Int_t firstBin, Int_t lastBin, Option_t* option) const
1958{
1959 /// method for performing projection
1960
1961 const char* expectedName = "_pr";
1962
1963 TString opt = option;
1964 TString cut;
1965 Int_t i1 = opt.Index("[");
1966 if(i1 >= 0) {
1967 Int_t i2 = opt.Index("]");
1968 cut = opt(i1, i2 - i1 + 1);
1969 }
1970 opt.ToLower(); // must be called after having parsed the cut name
1971 bool originalRange = opt.Contains("o");
1972
1973 Int_t firstXBin = fXaxis.GetFirst();
1974 Int_t lastXBin = fXaxis.GetLast();
1975
1976 if(firstXBin == 0 && lastXBin == 0) {
1977 firstXBin = 1;
1978 lastXBin = fXaxis.GetNbins();
1979 }
1980
1981 if(lastBin < firstBin && fYaxis.TestBit(TAxis::kAxisRange)) {
1982 firstBin = fYaxis.GetFirst();
1983 lastBin = fYaxis.GetLast();
1984 // For special case of TAxis::SetRange, when first == 1 and last
1985 // = N and the range bit has been set, the TAxis will return 0
1986 // for both.
1987 if(firstBin == 0 && lastBin == 0) {
1988 firstBin = 1;
1989 lastBin = fYaxis.GetNbins();
1990 }
1991 }
1992 if(firstBin < 0) {
1993 firstBin = 0;
1994 }
1995 if(lastBin < 0) {
1996 lastBin = fYaxis.GetLast() + 1;
1997 }
1998 if(lastBin > fYaxis.GetLast() + 1) {
1999 lastBin = fYaxis.GetLast() + 1;
2000 }
2001
2002 // Create the projection histogram
2003 char* pname = const_cast<char*>(name); // NOLINT(cppcoreguidelines-pro-type-const-cast)
2004 if(name != nullptr && strcmp(name, expectedName) == 0) {
2005 auto nch = strlen(GetName()) + 4;
2006 pname = new char[nch];
2007 snprintf(pname, nch, "%s%s", GetName(), name);
2008 }
2009 TH1D* h1 = nullptr;
2010 // check if histogram with identical name exist
2011 // if compatible reset and re-use previous histogram
2012 // (see https://savannah.cern.ch/bugs/?54340)
2013 TObject* h1obj = gROOT->FindObject(pname);
2014 if((h1obj != nullptr) && h1obj->InheritsFrom(TH1::Class())) {
2015 if(h1obj->IsA() != TH1D::Class()) {
2016 Error("DoProjection", "Histogram with name %s must be a TH1D and is a %s", name, h1obj->ClassName());
2017 return nullptr;
2018 }
2019 h1 = static_cast<TH1D*>(h1obj);
2020 // reset the existing histogram and set always the new binning for the axis
2021 // This avoid problems when the histogram already exists and the histograms is rebinned or its range has changed
2022 // (see https://savannah.cern.ch/bugs/?94101 or https://savannah.cern.ch/bugs/?95808 )
2023 h1->Reset();
2024 const TArrayD* xbins = fXaxis.GetXbins();
2025 if(xbins->fN == 0) {
2026 if(originalRange) {
2027 h1->SetBins(fXaxis.GetNbins(), fXaxis.GetXmin(), fXaxis.GetXmax());
2028 } else {
2029 h1->SetBins(lastXBin - firstXBin + 1, fXaxis.GetBinLowEdge(firstXBin), fXaxis.GetBinUpEdge(lastXBin));
2030 }
2031 } else {
2032 // case variable bins
2033 if(originalRange) {
2034 h1->SetBins(fXaxis.GetNbins(), xbins->fArray);
2035 } else {
2036 h1->SetBins(lastXBin - firstXBin + 1, &(xbins->fArray[firstXBin - 1]));
2037 }
2038 }
2039 }
2040 Int_t ncuts = 0;
2041 if(opt.Contains("[")) {
2042 const_cast<GHSym*>(this)->GetPainter(); // NOLINT(cppcoreguidelines-pro-type-const-cast)
2043 if(fPainter != nullptr) {
2044 ncuts = fPainter->MakeCuts(const_cast<char*>(cut.Data())); // NOLINT(cppcoreguidelines-pro-type-const-cast)
2045 }
2046 }
2047
2048 if(h1 == nullptr) {
2049 const TArrayD* bins = fXaxis.GetXbins();
2050 if(bins->fN == 0) {
2051 if(originalRange) {
2052 h1 = new TH1D(pname, GetTitle(), fXaxis.GetNbins(), fXaxis.GetXmin(), fXaxis.GetXmax());
2053 } else {
2054 h1 = new TH1D(pname, GetTitle(), lastXBin - firstXBin + 1, fXaxis.GetBinLowEdge(firstXBin),
2055 fXaxis.GetBinUpEdge(lastXBin));
2056 }
2057 } else {
2058 // case variable bins
2059 if(originalRange) {
2060 h1 = new TH1D(pname, GetTitle(), fXaxis.GetNbins(), bins->fArray);
2061 } else {
2062 h1 = new TH1D(pname, GetTitle(), lastXBin - firstXBin + 1, &(bins->fArray[firstXBin - 1]));
2063 }
2064 }
2065 if(opt.Contains("e") || (GetSumw2N() != 0)) {
2066 h1->Sumw2();
2067 }
2068 }
2069 if(pname != name) {
2070 delete[] pname;
2071 }
2072
2073 // Copy the axis attributes and the axis labels if needed.
2074 h1->GetXaxis()->ImportAttributes(&fXaxis);
2075 THashList* labels = const_cast<TAxis*>(&fXaxis)->GetLabels(); // NOLINT(cppcoreguidelines-pro-type-const-cast)
2076 if(labels != nullptr) {
2077 TIter iL(labels);
2078 TObjString* lb = nullptr;
2079 Int_t i = 1;
2080 while((lb = static_cast<TObjString*>(iL())) != nullptr) {
2081 h1->GetXaxis()->SetBinLabel(i, lb->String().Data());
2082 i++;
2083 }
2084 }
2085
2086 h1->SetLineColor(GetLineColor());
2087 h1->SetFillColor(GetFillColor());
2088 h1->SetMarkerColor(GetMarkerColor());
2089 h1->SetMarkerStyle(GetMarkerStyle());
2090
2091 // Fill the projected histogram
2092 Double_t totcont = 0;
2093 Bool_t computeErrors = h1->GetSumw2N() != 0;
2094
2095 // implement filling of projected histogram
2096 // xbin is bin number of xAxis (the projected axis). Loop is done on all bin of TH2 histograms
2097 // inbin is the axis being integrated. Loop is done only on the selected bins
2098 for(Int_t xbin = 0; xbin <= fXaxis.GetNbins() + 1; ++xbin) {
2099 Double_t err2 = 0.;
2100 Double_t cont = 0.;
2101 if(fXaxis.TestBit(TAxis::kAxisRange) && (xbin < firstXBin || xbin > lastXBin)) {
2102 continue;
2103 }
2104
2105 for(Int_t ybin = firstBin; ybin <= lastBin; ++ybin) {
2106 if(ncuts != 0) {
2107 if(!fPainter->IsInside(xbin, ybin)) {
2108 continue;
2109 }
2110 }
2111 // sum bin content and error if needed
2112 cont += GetCellContent(xbin, ybin);
2113 if(computeErrors) {
2114 Double_t exy = GetCellError(xbin, ybin);
2115 err2 += exy * exy;
2116 }
2117 }
2118 // find corresponding bin number in h1 for xbin
2119 Int_t binOut = h1->GetXaxis()->FindBin(fXaxis.GetBinCenter(xbin));
2120 h1->SetBinContent(binOut, cont);
2121 if(computeErrors) {
2122 h1->SetBinError(binOut, TMath::Sqrt(err2));
2123 }
2124 // sum all content
2125 totcont += cont;
2126 }
2127
2128 // check if we can re-use the original statistics from the previous histogram
2129 bool reuseStats = false;
2130 if((!fgStatOverflows && firstBin == 1 && lastBin == fYaxis.GetLast()) ||
2131 (fgStatOverflows && firstBin == 0 && lastBin == fYaxis.GetLast() + 1)) {
2132 reuseStats = true;
2133 } else {
2134 // also if total content match we can re-use
2135 double eps = 1.E-12;
2136 if(IsA() == GHSymF::Class()) {
2137 eps = 1.E-6;
2138 }
2139 if(fTsumw != 0 && TMath::Abs(fTsumw - totcont) < TMath::Abs(fTsumw) * eps) {
2140 reuseStats = true;
2141 }
2142 }
2143 if(ncuts != 0) {
2144 reuseStats = false;
2145 }
2146 // retrieve the statistics and set in projected histogram if we can re-use it
2147 bool reuseEntries = reuseStats;
2148 // can re-use entries if underflow/overflow are included
2149 reuseEntries &= static_cast<int>(firstBin == 0 && lastBin == fYaxis.GetLast() + 1);
2150 if(reuseStats) {
2151 std::array<Double_t, kNstat> stats;
2152 GetStats(stats.data());
2153 h1->PutStats(stats.data());
2154 } else {
2155 // the statistics is automatically recalulated since it is reset by the call to SetBinContent
2156 // we just need to set the entries since they have not been correctly calculated during the projection
2157 // we can only set them to the effective entries
2158 h1->SetEntries(h1->GetEffectiveEntries());
2159 }
2160 if(reuseEntries) {
2161 h1->SetEntries(fEntries);
2162 } else {
2163 // re-compute the entries
2164 // in case of error calculation (i.e. when Sumw2() is set)
2165 // use the effective entries for the entries
2166 // since this is the only way to estimate them
2167 Double_t entries = TMath::Floor(totcont + 0.5); // to avoid numerical rounding
2168 if(h1->GetSumw2N() != 0) {
2169 entries = h1->GetEffectiveEntries();
2170 }
2171 h1->SetEntries(entries);
2172 }
2173
2174 if(opt.Contains("d")) {
2175 TVirtualPad* padsav = gPad;
2176 TVirtualPad* pad = gROOT->GetSelectedPad();
2177 if(pad != nullptr) {
2178 pad->cd();
2179 }
2180 opt.Remove(opt.First("d"), 1);
2181 // remove also other options
2182 if(opt.Contains("e")) {
2183 opt.Remove(opt.First("e"), 1);
2184 }
2185 if(!gPad || !gPad->FindObject(h1)) {
2186 h1->Draw(opt);
2187 } else {
2188 h1->Paint(opt);
2189 }
2190 if(padsav != nullptr) {
2191 padsav->cd();
2192 }
2193 }
2194
2195 return h1;
2196}
2197
2198void GHSym::PutStats(Double_t* stats)
2199{
2200 // Replace current statistics with the values in array stats
2201 TH1::PutStats(stats);
2202 fTsumwy = stats[4];
2203 fTsumwy2 = stats[5];
2204 fTsumwxy = stats[6];
2205}
2206
2207GHSym* GHSym::Rebin2D(Int_t ngroup, const char* newname)
2208{
2209 // -*-*-*Rebin this histogram grouping ngroup/ngroup bins along the xaxis/yaxis together*-*-*-*-
2210 // =================================================================================
2211 // if newname is not blank a new temporary histogram hnew is created.
2212 // else the current histogram is modified (default)
2213 // The parameter ngroup indicates how many bins along the xaxis/yaxis of this
2214 // have to me merged into one bin of hnew
2215 // If the original histogram has errors stored (via Sumw2), the resulting
2216 // histograms has new errors correctly calculated.
2217 //
2218 // examples: if hpxpy is an existing GHSym histogram with 40 x 40 bins
2219 // hpxpy->Rebin2D(); // merges two bins along the xaxis and yaxis in one in hpxpy
2220 // // Carefull: previous contents of hpxpy are lost
2221 // hpxpy->Rebin2D(5); //merges five bins along the xaxisi and yaxis in one in hpxpy
2222 // GHSym* hnew = hpxpy->Rebin2D(5,"hnew"); // creates a new histogram hnew
2223 // // merging 5 bins of h1 along the xaxis and yaxis in one bin
2224 //
2225 // NOTE : If ngroup is not an exact divider of the number of bins,
2226 // along the xaxis/yaxis the top limit(s) of the rebinned histogram
2227 // is changed to the upper edge of the bin=newbins*ngroup
2228 // and the corresponding bins are added to
2229 // the overflow bin.
2230 // Statistics will be recomputed from the new bin contents.
2231
2232 Int_t nbins = fXaxis.GetNbins();
2233 Double_t min = fXaxis.GetXmin();
2234 Double_t max = fXaxis.GetXmax();
2235 if((ngroup <= 0) || (ngroup > nbins)) {
2236 Error("Rebin", "Illegal value of ngroup=%d", ngroup);
2237 return nullptr;
2238 }
2239
2240 Int_t newbins = nbins / ngroup;
2241
2242 // Save old bin contents into a new array
2243 Double_t entries = fEntries;
2244 auto* oldBins = new Double_t[(nbins + 2) * (nbins + 3) / 2];
2245 for(Int_t xbin = 0; xbin < nbins + 2; xbin++) {
2246 for(Int_t ybin = 0; ybin <= xbin; ybin++) {
2247 Int_t bin = GetBin(xbin, ybin);
2248 oldBins[bin] = GetBinContent(bin);
2249 }
2250 }
2251 Double_t* oldErrors = nullptr;
2252 if(fSumw2.fN != 0) {
2253 oldErrors = new Double_t[(nbins + 2) * (nbins + 3) / 2];
2254 for(Int_t xbin = 0; xbin < nbins + 2; xbin++) {
2255 for(Int_t ybin = 0; ybin <= xbin; ybin++) {
2256 Int_t bin = GetBin(xbin, ybin);
2257 oldErrors[bin] = GetBinError(bin);
2258 }
2259 }
2260 }
2261
2262 // create a clone of the old histogram if newname is specified
2263 GHSym* hnew = this;
2264 if((newname != nullptr) && (strlen(newname) != 0u)) {
2265 hnew = static_cast<GHSym*>(Clone());
2266 hnew->SetName(newname);
2267 }
2268
2269 // save original statistics
2270 std::array<Double_t, kNstat> stats;
2271 GetStats(stats.data());
2272 bool resetStat = false;
2273
2274 // change axis specs and rebuild bin contents array
2275 if(newbins * ngroup != nbins) {
2276 max = fXaxis.GetBinUpEdge(newbins * ngroup);
2277 resetStat = true; // stats must be reset because top bins will be moved to overflow bin
2278 }
2279 // save the TAttAxis members (reset by SetBins) for x axis
2280 Int_t nXdivisions = fXaxis.GetNdivisions();
2281 Color_t xAxisColor = fXaxis.GetAxisColor();
2282 Color_t xLabelColor = fXaxis.GetLabelColor();
2283 Style_t xLabelFont = fXaxis.GetLabelFont();
2284 Float_t xLabelOffset = fXaxis.GetLabelOffset();
2285 Float_t xLabelSize = fXaxis.GetLabelSize();
2286 Float_t xTickLength = fXaxis.GetTickLength();
2287 Float_t xTitleOffset = fXaxis.GetTitleOffset();
2288 Float_t xTitleSize = fXaxis.GetTitleSize();
2289 Color_t xTitleColor = fXaxis.GetTitleColor();
2290 Style_t xTitleFont = fXaxis.GetTitleFont();
2291 // save the TAttAxis members (reset by SetBins) for y axis
2292 Int_t nYdivisions = fYaxis.GetNdivisions();
2293 Color_t yAxisColor = fYaxis.GetAxisColor();
2294 Color_t yLabelColor = fYaxis.GetLabelColor();
2295 Style_t yLabelFont = fYaxis.GetLabelFont();
2296 Float_t yLabelOffset = fYaxis.GetLabelOffset();
2297 Float_t yLabelSize = fYaxis.GetLabelSize();
2298 Float_t yTickLength = fYaxis.GetTickLength();
2299 Float_t yTitleOffset = fYaxis.GetTitleOffset();
2300 Float_t yTitleSize = fYaxis.GetTitleSize();
2301 Color_t yTitleColor = fYaxis.GetTitleColor();
2302 Style_t yTitleFont = fYaxis.GetTitleFont();
2303
2304 // copy merged bin contents (ignore under/overflows)
2305 if(ngroup != 1) {
2306 if(fXaxis.GetXbins()->GetSize() > 0 || fYaxis.GetXbins()->GetSize() > 0) {
2307 // variable bin sizes in x or y, don't treat both cases separately
2308 auto* bins = new Double_t[newbins + 1];
2309 for(Int_t i = 0; i <= newbins; ++i) {
2310 bins[i] = fXaxis.GetBinLowEdge(1 + i * ngroup);
2311 }
2312 hnew->SetBins(newbins, bins, newbins, bins); // changes also errors array (if any)
2313 delete[] bins;
2314 } else {
2315 hnew->SetBins(newbins, min, max, newbins, min, max); // changes also errors array
2316 }
2317
2318 Int_t oldxbin = 1;
2319 Int_t oldybin = 1;
2320 for(Int_t xbin = 1; xbin <= newbins; ++xbin) {
2321 oldybin = 1;
2322 for(Int_t ybin = 1; ybin <= xbin; ++ybin) {
2323 Double_t binContent = 0.;
2324 Double_t binError = 0.;
2325 for(Int_t i = 0; i < ngroup; ++i) {
2326 if(oldxbin + i > nbins) {
2327 break;
2328 }
2329 for(Int_t j = 0; j < ngroup; ++j) {
2330 if(oldybin + j > nbins) {
2331 break;
2332 }
2333 // get global bin (same conventions as in GHSym::GetBin(xbin,ybin)
2334 Int_t bin = 0;
2335 if(oldybin + j <= oldxbin + i) {
2336 bin = oldxbin + i + (oldybin + j) * (2 * fXaxis.GetNbins() - (oldybin + j) + 3) / 2;
2337 } else {
2338 bin = oldybin + j + (oldxbin + i) * (2 * fXaxis.GetNbins() - (oldxbin + i) + 3) / 2;
2339 }
2340 binContent += oldBins[bin];
2341 if(oldErrors != nullptr) {
2342 binError += oldErrors[bin] * oldErrors[bin];
2343 }
2344 }
2345 }
2346 hnew->SetBinContent(xbin, ybin, binContent);
2347 if(oldErrors != nullptr) {
2348 hnew->SetBinError(xbin, ybin, TMath::Sqrt(binError));
2349 }
2350 oldybin += ngroup;
2351 }
2352 oldxbin += ngroup;
2353 }
2354
2355 // Recompute correct underflows and overflows.
2356
2357 // copy old underflow bin in x and y (0,0)
2358 hnew->SetBinContent(0, 0, oldBins[0]);
2359 if(oldErrors != nullptr) {
2360 hnew->SetBinError(0, 0, oldErrors[0]);
2361 }
2362
2363 // calculate new overflow bin in x and y (newbins+1,newbins+1)
2364 Double_t binContent = 0.;
2365 Double_t binError = 0.;
2366 for(Int_t xbin = oldxbin; xbin <= nbins + 1; ++xbin) {
2367 for(Int_t ybin = oldybin; ybin <= xbin; ++ybin) {
2368 Int_t bin = xbin + ybin * (2 * nbins - ybin + 3) / 2;
2369 binContent += oldBins[bin];
2370 if(oldErrors != nullptr) {
2371 binError += oldErrors[bin] * oldErrors[bin];
2372 }
2373 }
2374 }
2375 hnew->SetBinContent(newbins + 1, newbins + 1, binContent);
2376 if(oldErrors != nullptr) {
2377 hnew->SetBinError(newbins + 1, newbins + 1, TMath::Sqrt(binError));
2378 }
2379
2380 // calculate new underflow bin in x and overflow in y (0,newbins+1)
2381 binContent = 0.;
2382 binError = 0.;
2383 for(Int_t ybin = oldybin; ybin <= nbins + 1; ++ybin) {
2384 Int_t bin = ybin * (2 * nbins - ybin + 3) / 2;
2385 binContent += oldBins[bin];
2386 if(oldErrors != nullptr) {
2387 binError += oldErrors[bin] * oldErrors[bin];
2388 }
2389 }
2390 hnew->SetBinContent(0, newbins + 1, binContent);
2391 if(oldErrors != nullptr) {
2392 hnew->SetBinError(0, newbins + 1, TMath::Sqrt(binError));
2393 }
2394
2395 // calculate new overflow bin in x and underflow in y (newbins+1,0)
2396 binContent = 0.;
2397 binError = 0.;
2398 for(Int_t xbin = oldxbin; xbin <= nbins + 1; ++xbin) {
2399 Int_t bin = xbin;
2400 binContent += oldBins[bin];
2401 if(oldErrors != nullptr) {
2402 binError += oldErrors[bin] * oldErrors[bin];
2403 }
2404 }
2405 hnew->SetBinContent(newbins + 1, 0, binContent);
2406 if(oldErrors != nullptr) {
2407 hnew->SetBinError(newbins + 1, 0, TMath::Sqrt(binError));
2408 }
2409
2410 // recompute under/overflow contents in y for the new x bins
2411 Int_t oldxbin2 = 1;
2412 for(Int_t xbin = 1; xbin <= newbins; ++xbin) {
2413 Double_t binContent0 = 0.;
2414 Double_t binContent2 = 0.;
2415 Double_t binError0 = 0.;
2416 Double_t binError2 = 0.;
2417 for(Int_t i = 0; i < ngroup; ++i) {
2418 if(oldxbin2 + i > nbins) {
2419 break;
2420 }
2421 // old underflow bin (in y)
2422 Int_t ufbin = oldxbin2 + i;
2423 binContent0 += oldBins[ufbin];
2424 if(oldErrors != nullptr) {
2425 binError0 += oldErrors[ufbin] * oldErrors[ufbin];
2426 }
2427 for(Int_t ybin = oldybin; ybin <= nbins + 1; ++ybin) {
2428 // old overflow bin (in y)
2429 Int_t ofbin = ufbin + ybin * (nbins + 2);
2430 binContent2 += oldBins[ofbin];
2431 if(oldErrors != nullptr) {
2432 binError2 += oldErrors[ofbin] * oldErrors[ofbin];
2433 }
2434 }
2435 }
2436 hnew->SetBinContent(xbin, 0, binContent0);
2437 hnew->SetBinContent(xbin, newbins + 1, binContent2);
2438 if(oldErrors != nullptr) {
2439 hnew->SetBinError(xbin, 0, TMath::Sqrt(binError0));
2440 hnew->SetBinError(xbin, newbins + 1, TMath::Sqrt(binError2));
2441 }
2442 oldxbin2 += ngroup;
2443 }
2444
2445 // recompute under/overflow contents in x for the new y bins
2446 Int_t oldybin2 = 1;
2447 for(Int_t ybin = 1; ybin <= newbins; ++ybin) {
2448 Double_t binContent0 = 0.;
2449 Double_t binContent2 = 0.;
2450 Double_t binError0 = 0.;
2451 Double_t binError2 = 0.;
2452 for(Int_t i = 0; i < ngroup; ++i) {
2453 if(oldybin2 + i > nbins) {
2454 break;
2455 }
2456 // old underflow bin (in x)
2457 Int_t ufbin = (oldybin2 + i) * (nbins + 2);
2458 binContent0 += oldBins[ufbin];
2459 if(oldErrors != nullptr) {
2460 binError0 += oldErrors[ufbin] * oldErrors[ufbin];
2461 }
2462 for(Int_t xbin = oldxbin; xbin <= nbins + 1; ++xbin) {
2463 Int_t ofbin = ufbin + xbin;
2464 binContent2 += oldBins[ofbin];
2465 if(oldErrors != nullptr) {
2466 binError2 += oldErrors[ofbin] * oldErrors[ofbin];
2467 }
2468 }
2469 }
2470 hnew->SetBinContent(0, ybin, binContent0);
2471 hnew->SetBinContent(newbins + 1, ybin, binContent2);
2472 if(oldErrors != nullptr) {
2473 hnew->SetBinError(0, ybin, TMath::Sqrt(binError0));
2474 hnew->SetBinError(newbins + 1, ybin, TMath::Sqrt(binError2));
2475 }
2476 oldybin2 += ngroup;
2477 }
2478 }
2479
2480 // Restore x axis attributes
2481 fXaxis.SetNdivisions(nXdivisions);
2482 fXaxis.SetAxisColor(xAxisColor);
2483 fXaxis.SetLabelColor(xLabelColor);
2484 fXaxis.SetLabelFont(xLabelFont);
2485 fXaxis.SetLabelOffset(xLabelOffset);
2486 fXaxis.SetLabelSize(xLabelSize);
2487 fXaxis.SetTickLength(xTickLength);
2488 fXaxis.SetTitleOffset(xTitleOffset);
2489 fXaxis.SetTitleSize(xTitleSize);
2490 fXaxis.SetTitleColor(xTitleColor);
2491 fXaxis.SetTitleFont(xTitleFont);
2492 // Restore y axis attributes
2493 fYaxis.SetNdivisions(nYdivisions);
2494 fYaxis.SetAxisColor(yAxisColor);
2495 fYaxis.SetLabelColor(yLabelColor);
2496 fYaxis.SetLabelFont(yLabelFont);
2497 fYaxis.SetLabelOffset(yLabelOffset);
2498 fYaxis.SetLabelSize(yLabelSize);
2499 fYaxis.SetTickLength(yTickLength);
2500 fYaxis.SetTitleOffset(yTitleOffset);
2501 fYaxis.SetTitleSize(yTitleSize);
2502 fYaxis.SetTitleColor(yTitleColor);
2503 fYaxis.SetTitleFont(yTitleFont);
2504
2505 // restore statistics and entries modified by SetBinContent
2506 hnew->SetEntries(entries);
2507 if(!resetStat) {
2508 hnew->PutStats(stats.data());
2509 }
2510
2511 delete[] oldBins;
2512 delete[] oldErrors;
2513 return hnew;
2514}
2515
2516void GHSym::Reset(Option_t* option)
2517{
2518 //*-*-*-*-*-*-*-*Reset this histogram: contents, errors, etc*-*-*-*-*-*-*-*
2519 //*-* ===========================================
2520
2521 TH1::Reset(option);
2522 TString opt = option;
2523 opt.ToUpper();
2524
2525 if(opt.Contains("ICE") && !opt.Contains("S")) {
2526 return;
2527 }
2528 fTsumwy = 0;
2529 fTsumwy2 = 0;
2530 fTsumwxy = 0;
2531 if(fMatrix != nullptr) {
2532 try {
2533 delete fMatrix;
2534 fMatrix = nullptr;
2535 } catch(std::exception& e) {
2536 fMatrix = nullptr;
2537 }
2538 }
2539}
2540
2541void GHSym::SetCellContent(Int_t binx, Int_t biny, Double_t content)
2542{
2543 SetBinContent(GetBin(binx, biny), content);
2544}
2545
2546void GHSym::SetCellError(Int_t binx, Int_t biny, Double_t content)
2547{
2548 SetBinError(GetBin(binx, biny), content);
2549}
2550
2552{
2553 // When the mouse is moved in a pad containing a 2-d view of this histogram
2554 // a second canvas shows the projection along X corresponding to the
2555 // mouse position along Y.
2556 // To stop the generation of the projections, delete the canvas
2557 // containing the projection.
2558
2559 GetPainter();
2560 if(fPainter != nullptr) {
2561 fPainter->SetShowProjection("x", nbins);
2562 }
2563}
2564
2566{
2567 // When the mouse is moved in a pad containing a 2-d view of this histogram
2568 // a second canvas shows the projection along Y corresponding to the
2569 // mouse position along X.
2570 // To stop the generation of the projections, delete the canvas
2571 // containing the projection.
2572
2573 GetPainter();
2574 if(fPainter != nullptr) {
2575 fPainter->SetShowProjection("y", nbins);
2576 }
2577}
2578
2579TH1* GHSym::ShowBackground(Int_t niter, Option_t* option)
2580{
2581 // This function calculates the background spectrum in this histogram.
2582 // The background is returned as a histogram.
2583 // to be implemented (may be)
2584
2585 return reinterpret_cast<TH1*>(gROOT->ProcessLineFast( // NOLINT(performance-no-int-to-ptr)
2586 Form(R"(TSpectrum2::StaticBackground((TH1*)0x%lx,%d,"%s"))", reinterpret_cast<ULong_t>(this), niter, option)));
2587}
2588
2589Int_t GHSym::ShowPeaks(Double_t sigma, Option_t* option, Double_t threshold)
2590{
2591 // Interface to TSpectrum2::Search
2592 // the function finds peaks in this histogram where the width is > sigma
2593 // and the peak maximum greater than threshold*maximum bin content of this.
2594 // for more detauils see TSpectrum::Search.
2595 // note the difference in the default value for option compared to TSpectrum2::Search
2596 // option="" by default (instead of "goff")
2597
2598 return static_cast<Int_t>(gROOT->ProcessLineFast(
2599 Form(R"(TSpectrum2::StaticSearch((TH1*)0x%lx,%g,"%s",%g))", reinterpret_cast<ULong_t>(this), sigma, option, threshold)));
2600}
2601
2602void GHSym::Smooth(Int_t ntimes, Option_t* option)
2603{
2604 // Smooth bin contents of this 2-d histogram using kernel algorithms
2605 // similar to the ones used in the raster graphics community.
2606 // Bin contents in the active range are replaced by their smooth values.
2607 // If Errors are defined via Sumw2, they are also scaled and computed.
2608 // However, note the resulting errors will be correlated between different-bins, so
2609 // the errors should not be used blindly to perform any calculation involving several bins,
2610 // like fitting the histogram. One would need to compute also the bin by bin correlation matrix.
2611 //
2612 // 3 kernels are proposed k5a, k5b and k3a.
2613 // k5a and k5b act on 5x5 cells (i-2,i-1,i,i+1,i+2, and same for j)
2614 // k5b is a bit more stronger in smoothing
2615 // k3a acts only on 3x3 cells (i-1,i,i+1, and same for j).
2616 // By default the kernel "k5a" is used. You can select the kernels "k5b" or "k3a"
2617 // via the option argument.
2618 // If TAxis::SetRange has been called on the x or/and y axis, only the bins
2619 // in the specified range are smoothed.
2620 // In the current implementation if the first argument is not used (default value=1).
2621 //
2622 // implementation by David McKee (dmckee@bama.ua.edu). Extended by Rene Brun
2623
2624 std::array<std::array<Double_t, 5>, 5> k5a = {{{0, 0, 1, 0, 0}, {0, 2, 2, 2, 0}, {1, 2, 5, 2, 1}, {0, 2, 2, 2, 0}, {0, 0, 1, 0, 0}}};
2625 std::array<std::array<Double_t, 5>, 5> k5b = {{{0, 1, 2, 1, 0}, {1, 2, 4, 2, 1}, {2, 4, 8, 4, 2}, {1, 2, 4, 2, 1}, {0, 1, 2, 1, 0}}};
2626 std::array<std::array<Double_t, 3>, 3> k3a = {{{0, 1, 0}, {1, 2, 1}, {0, 1, 0}}};
2627
2628 if(ntimes > 1) {
2629 Warning("Smooth", "Currently only ntimes=1 is supported");
2630 }
2631 TString opt = option;
2632 opt.ToLower();
2633 Int_t ksize_x = k5a.size();
2634 Int_t ksize_y = k5a.size();
2635 Double_t* kernel = k5a.data()->data();
2636 if(opt.Contains("k5b")) {
2637 kernel = k5b.data()->data();
2638 }
2639 if(opt.Contains("k3a")) {
2640 kernel = k3a.data()->data();
2641 ksize_x = k3a.size();
2642 ksize_y = k3a.size();
2643 }
2644
2645 // find i,j ranges
2646 Int_t ifirst = fXaxis.GetFirst();
2647 Int_t ilast = fXaxis.GetLast();
2648 Int_t jfirst = fYaxis.GetFirst();
2649 Int_t jlast = fYaxis.GetLast();
2650
2651 // Determine the size of the bin buffer(s) needed
2652 Double_t nentries = fEntries;
2653 Int_t nx = GetNbinsX();
2654 Int_t ny = GetNbinsY();
2655 Int_t bufSize = (nx + 2) * (ny + 2);
2656 auto* buf = new Double_t[bufSize];
2657 Double_t* ebuf = nullptr;
2658 if(fSumw2.fN != 0) {
2659 ebuf = new Double_t[bufSize];
2660 }
2661
2662 // Copy all the data to the temporary buffers
2663 for(Int_t i = ifirst; i <= ilast; ++i) {
2664 for(Int_t j = jfirst; j <= jlast; ++j) {
2665 Int_t bin = GetBin(i, j);
2666 buf[bin] = GetBinContent(bin);
2667 if(ebuf != nullptr) {
2668 ebuf[bin] = GetBinError(bin);
2669 }
2670 }
2671 }
2672
2673 // Kernel tail sizes (kernel sizes must be odd for this to work!)
2674 Int_t x_push = (ksize_x - 1) / 2;
2675 Int_t y_push = (ksize_y - 1) / 2;
2676
2677 // main work loop
2678 for(Int_t i = ifirst; i <= ilast; ++i) {
2679 for(Int_t j = jfirst; j <= jlast; ++j) {
2680 Double_t content = 0.0;
2681 Double_t error = 0.0;
2682 Double_t norm = 0.0;
2683
2684 for(Int_t n = 0; n < ksize_x; ++n) {
2685 for(Int_t m = 0; m < ksize_y; ++m) {
2686 Int_t xb = i + (n - x_push);
2687 Int_t yb = j + (m - y_push);
2688 if((xb >= 1) && (xb <= nx) && (yb >= 1) && (yb <= ny)) {
2689 Int_t bin = GetBin(xb, yb);
2690 Double_t k = kernel[n * ksize_y + m];
2691 if(k != 0.0) {
2692 norm += k;
2693 content += k * buf[bin];
2694 if(ebuf != nullptr) {
2695 error += k * k * ebuf[bin] * ebuf[bin];
2696 }
2697 }
2698 }
2699 }
2700 }
2701
2702 if(norm != 0.0) {
2703 SetBinContent(i, j, content / norm);
2704 if(ebuf != nullptr) {
2705 error /= (norm * norm);
2706 SetBinError(i, j, sqrt(error));
2707 }
2708 }
2709 }
2710 }
2711 fEntries = nentries;
2712
2713 delete[] buf;
2714 delete[] ebuf;
2715}
2716
2717//------------------------------------------------------------
2718// GHSymF methods (float = four bytes per cell)
2719//------------------------------------------------------------
2720
2722{
2723 SetBinsLength(9);
2724 if(fgDefaultSumw2) {
2725 Sumw2();
2726 }
2727}
2728
2729GHSymF::GHSymF(const char* name, const char* title, Int_t nbins, Double_t low, Double_t up)
2730 : GHSym(name, title, nbins, low, up)
2731{
2732 TArrayF::Set(fNcells);
2733 if(fgDefaultSumw2) {
2734 Sumw2();
2735 }
2736
2737 if(low >= up) {
2738 SetBuffer(fgBufferSize);
2739 }
2740}
2741
2742GHSymF::GHSymF(const char* name, const char* title, Int_t nbins, const Double_t* bins) : GHSym(name, title, nbins, bins)
2743{
2744 TArrayF::Set(fNcells);
2745 if(fgDefaultSumw2) {
2746 Sumw2();
2747 }
2748}
2749
2750GHSymF::GHSymF(const char* name, const char* title, Int_t nbins, const Float_t* bins) : GHSym(name, title, nbins, bins)
2751{
2752 TArrayF::Set(fNcells);
2753 if(fgDefaultSumw2) {
2754 Sumw2();
2755 }
2756}
2757
2759 : GHSym(rhs), TArrayF(rhs)
2760{
2761 rhs.Copy(*this);
2762}
2763
2764GHSymF::GHSymF(GHSymF&& rhs) noexcept
2765 : GHSym(std::move(rhs)), TArrayF(rhs)
2766{
2767 rhs.Copy(*this);
2768}
2769
2770GHSymF::~GHSymF() = default;
2771
2772TH2F* GHSymF::GetMatrix(bool force)
2773{
2774 if(Matrix() != nullptr && !force) {
2775 return static_cast<TH2F*>(Matrix());
2776 }
2777 if(force) {
2778 delete Matrix();
2779 }
2780
2781 Matrix(new TH2F(Form("%s_mat", GetName()), GetTitle(), fXaxis.GetNbins(), fXaxis.GetXmin(), fXaxis.GetXmax(),
2782 fYaxis.GetNbins(), fYaxis.GetXmin(), fYaxis.GetXmax()));
2783 // copy cell contents (including all overflow and underflow cells)
2784 for(int i = 0; i < fXaxis.GetNbins() + 2; ++i) {
2785 for(int j = 0; j < fXaxis.GetNbins() + 2; ++j) {
2786 Matrix()->SetBinContent(i, j, GetBinContent(i, j));
2787 }
2788 }
2789 return static_cast<TH2F*>(Matrix());
2790}
2791
2792void GHSymF::Copy(TObject& rhs) const
2793{
2794 GHSym::Copy(static_cast<GHSymF&>(rhs));
2795}
2796
2797TH1* GHSymF::DrawCopy(Option_t* option, const char* name_postfix) const
2798{
2799 // Draw copy.
2800
2801 TString opt = option;
2802 opt.ToLower();
2803 if(gPad != nullptr && !opt.Contains("same")) {
2804 gPad->Clear();
2805 }
2806 TString newName = (name_postfix) != nullptr ? TString::Format("%s%s", GetName(), name_postfix) : "";
2807 TH1* newth1 = static_cast<TH1*>(Clone(newName));
2808 newth1->SetDirectory(nullptr);
2809 newth1->SetBit(kCanDelete);
2810 newth1->AppendPad(option);
2811 return newth1;
2812}
2813
2814Double_t GHSymF::GetBinContent(Int_t bin) const
2815{
2816 // Get bin content.
2817
2818 if(fBuffer != nullptr) {
2819 const_cast<GHSymF*>(this)->BufferEmpty(); // NOLINT(cppcoreguidelines-pro-type-const-cast)
2820 }
2821 if(bin < 0) {
2822 bin = 0;
2823 }
2824 if(bin >= fNcells) {
2825 bin = fNcells - 1;
2826 }
2827 if(fArray == nullptr) {
2828 return 0;
2829 }
2830 return static_cast<Double_t>(fArray[bin]);
2831}
2832
2833void GHSymF::Reset(Option_t* option)
2834{
2835 //*-*-*-*-*-*-*-*Reset this histogram: contents, errors, etc*-*-*-*-*-*-*-*
2836 //*-* ===========================================
2837
2838 GHSym::Reset(option);
2839 TArrayF::Reset();
2840}
2841
2842void GHSymF::SetBinContent(Int_t bin, Double_t content)
2843{
2844 // Set bin content
2845 fEntries++;
2846 fTsumw = 0;
2847 if(bin < 0) {
2848 return;
2849 }
2850 if(bin >= fNcells) {
2851 return;
2852 }
2853 fArray[bin] = static_cast<Float_t>(content);
2854}
2855
2857{
2858 // Set total number of bins including under/overflow
2859 // Reallocate bin contents array
2860
2861 if(n < 0) {
2862 n = (fXaxis.GetNbins() + 2) * (fYaxis.GetNbins() + 2);
2863 }
2864 fNcells = n;
2865 TArrayF::Set(n);
2866}
2867
2869{
2870 // Operator =
2871
2872 if(this != &h1) {
2873 h1.Copy(*this);
2874 }
2875 return *this;
2876}
2877
2879{
2880 // Operator =
2881
2882 if(this != &h1) {
2883 h1.Copy(*this);
2884 }
2885 return *this;
2886}
2887
2888GHSymF operator*(Float_t c1, GHSymF& h1)
2889{
2890 // Operator *
2891
2892 GHSymF hnew = h1;
2893 hnew.Scale(c1);
2894 hnew.SetDirectory(nullptr);
2895 return hnew;
2896}
2897
2899{
2900 // Operator +
2901
2902 GHSymF hnew = h1;
2903 hnew.Add(&h2, 1);
2904 hnew.SetDirectory(nullptr);
2905 return hnew;
2906}
2907
2909{
2910 // Operator -
2911
2912 GHSymF hnew = h1;
2913 hnew.Add(&h2, -1);
2914 hnew.SetDirectory(nullptr);
2915 return hnew;
2916}
2917
2919{
2920 // Operator *
2921
2922 GHSymF hnew = h1;
2923 hnew.Multiply(&h2);
2924 hnew.SetDirectory(nullptr);
2925 return hnew;
2926}
2927
2929{
2930 // Operator /
2931
2932 GHSymF hnew = h1;
2933 hnew.Divide(&h2);
2934 hnew.SetDirectory(nullptr);
2935 return hnew;
2936}
2937
2938//------------------------------------------------------------
2939// GHSymD methods (double = eight bytes per cell)
2940//------------------------------------------------------------
2941
2943{
2944 SetBinsLength(9);
2945 if(fgDefaultSumw2) {
2946 Sumw2();
2947 }
2948}
2949
2950GHSymD::GHSymD(const char* name, const char* title, Int_t nbins, Double_t low, Double_t up)
2951 : GHSym(name, title, nbins, low, up)
2952{
2953 TArrayD::Set(fNcells);
2954 if(fgDefaultSumw2) {
2955 Sumw2();
2956 }
2957
2958 if(low >= up) {
2959 SetBuffer(fgBufferSize);
2960 }
2961}
2962
2963GHSymD::GHSymD(const char* name, const char* title, Int_t nbins, const Double_t* bins) : GHSym(name, title, nbins, bins)
2964{
2965 TArrayD::Set(fNcells);
2966 if(fgDefaultSumw2) {
2967 Sumw2();
2968 }
2969}
2970
2971GHSymD::GHSymD(const char* name, const char* title, Int_t nbins, const Float_t* bins) : GHSym(name, title, nbins, bins)
2972{
2973 TArrayD::Set(fNcells);
2974 if(fgDefaultSumw2) {
2975 Sumw2();
2976 }
2977}
2978
2980 : GHSym(rhs), TArrayD(rhs)
2981{
2982 rhs.Copy(*this);
2983}
2984
2985GHSymD::GHSymD(GHSymD&& rhs) noexcept
2986 : GHSym(std::move(rhs)), TArrayD(rhs)
2987{
2988 rhs.Copy(*this);
2989}
2990
2991GHSymD::~GHSymD() = default;
2992
2993TH2D* GHSymD::GetMatrix(bool force)
2994{
2995 if(Matrix() != nullptr && !force) {
2996 return static_cast<TH2D*>(Matrix());
2997 }
2998 if(force) {
2999 delete Matrix();
3000 }
3001
3002 Matrix(new TH2D(Form("%s_mat", GetName()), Form("%s;%s;%s", GetTitle(), fXaxis.GetTitle(), fYaxis.GetTitle()),
3003 fXaxis.GetNbins(), fXaxis.GetXmin(), fXaxis.GetXmax(), fYaxis.GetNbins(), fYaxis.GetXmin(),
3004 fYaxis.GetXmax()));
3005 // copy cell contents (including all overflow and underflow cells)
3006 for(int i = 0; i < fXaxis.GetNbins() + 2; ++i) {
3007 for(int j = 0; j < fXaxis.GetNbins() + 2; ++j) {
3008 Matrix()->SetBinContent(i, j, GetBinContent(i, j));
3009 }
3010 }
3011 return static_cast<TH2D*>(Matrix());
3012}
3013
3014void GHSymD::Copy(TObject& rhs) const
3015{
3016 GHSym::Copy(static_cast<GHSymD&>(rhs));
3017}
3018
3019TH1* GHSymD::DrawCopy(Option_t* option, const char* name_postfix) const
3020{
3021 // Draw copy.
3022
3023 TString opt = option;
3024 opt.ToLower();
3025 if(gPad != nullptr && !opt.Contains("same")) {
3026 gPad->Clear();
3027 }
3028 TString newName = (name_postfix) != nullptr ? TString::Format("%s%s", GetName(), name_postfix) : "";
3029 TH1* newth1 = static_cast<TH1*>(Clone(newName));
3030 newth1->SetDirectory(nullptr);
3031 newth1->SetBit(kCanDelete);
3032 newth1->AppendPad(option);
3033 return newth1;
3034}
3035
3036Double_t GHSymD::GetBinContent(Int_t bin) const
3037{
3038 // Get bin content.
3039 if(fBuffer != nullptr) {
3040 const_cast<GHSymD*>(this)->BufferEmpty(); // NOLINT(cppcoreguidelines-pro-type-const-cast)
3041 }
3042 if(bin < 0) {
3043 bin = 0;
3044 }
3045 if(bin >= fNcells) {
3046 bin = fNcells - 1;
3047 }
3048 if(fArray == nullptr) {
3049 return 0;
3050 }
3051 return static_cast<Double_t>(fArray[bin]);
3052}
3053
3054void GHSymD::Reset(Option_t* option)
3055{
3056 //*-*-*-*-*-*-*-*Reset this histogram: contents, errors, etc*-*-*-*-*-*-*-*
3057 //*-* ===========================================
3058
3059 GHSym::Reset(option);
3060 TArrayD::Reset();
3061}
3062
3063void GHSymD::SetBinContent(Int_t bin, Double_t content)
3064{
3065 // Set bin content
3066 fEntries++;
3067 fTsumw = 0;
3068 if(bin < 0) {
3069 return;
3070 }
3071 if(bin >= fNcells) {
3072 return;
3073 }
3074 fArray[bin] = static_cast<Float_t>(content);
3075}
3076
3078{
3079 // Set total number of bins including under/overflow
3080 // Reallocate bin contents array
3081
3082 if(n < 0) {
3083 n = (fXaxis.GetNbins() + 2) * (fYaxis.GetNbins() + 2);
3084 }
3085 fNcells = n;
3086 TArrayD::Set(n);
3087}
3088
3090{
3091 // Operator =
3092
3093 if(this != &h1) {
3094 h1.Copy(*this);
3095 }
3096 return *this;
3097}
3098
3100{
3101 // Operator =
3102
3103 if(this != &h1) {
3104 h1.Copy(*this);
3105 }
3106 return *this;
3107}
3108
3109GHSymD operator*(Float_t c1, GHSymD& h1)
3110{
3111 // Operator *
3112
3113 GHSymD hnew = h1;
3114 hnew.Scale(c1);
3115 hnew.SetDirectory(nullptr);
3116 return hnew;
3117}
3118
3120{
3121 // Operator +
3122
3123 GHSymD hnew = h1;
3124 hnew.Add(&h2, 1);
3125 hnew.SetDirectory(nullptr);
3126 return hnew;
3127}
3128
3130{
3131 // Operator -
3132
3133 GHSymD hnew = h1;
3134 hnew.Add(&h2, -1);
3135 hnew.SetDirectory(nullptr);
3136 return hnew;
3137}
3138
3140{
3141 // Operator *
3142
3143 GHSymD hnew = h1;
3144 hnew.Multiply(&h2);
3145 hnew.SetDirectory(nullptr);
3146 return hnew;
3147}
3148
3150{
3151 // Operator /
3152
3153 GHSymD hnew = h1;
3154 hnew.Divide(&h2);
3155 hnew.SetDirectory(nullptr);
3156 return hnew;
3157}
3158
3159//------------------------------------------------------------
GHSymF operator/(GHSymF &h1, GHSymF &h2)
Definition GHSym.cxx:2928
GHSymF operator+(GHSymF &h1, GHSymF &h2)
Definition GHSym.cxx:2898
GHSymF operator-(GHSymF &h1, GHSymF &h2)
Definition GHSym.cxx:2908
GHSymF operator*(Float_t c1, GHSymF &h1)
Definition GHSym.cxx:2888
TH1D * hist
Definition UserFillObj.h:3
Double_t GetBinContent(Int_t bin) const override
Definition GHSym.cxx:3036
TH1 * DrawCopy(Option_t *option="", const char *name_postfix="_copy") const override
Definition GHSym.cxx:3019
GHSymD()
Definition GHSym.cxx:2942
TH2D * GetMatrix(bool force=false)
Definition GHSym.cxx:2993
void Reset(Option_t *option="") override
Definition GHSym.cxx:3054
void Copy(TObject &rh) const override
Definition GHSym.cxx:3014
void SetBinsLength(Int_t n=-1) override
Definition GHSym.cxx:3077
GHSymD & operator=(const GHSymD &h1)
Definition GHSym.cxx:3089
void SetBinContent(Int_t bin, Double_t content) override
Definition GHSym.cxx:3063
Double_t GetBinContent(Int_t bin) const override
Definition GHSym.cxx:2814
void SetBinsLength(Int_t n=-1) override
Definition GHSym.cxx:2856
GHSymF & operator=(const GHSymF &h1)
Definition GHSym.cxx:2868
void Reset(Option_t *option="") override
Definition GHSym.cxx:2833
void SetBinContent(Int_t bin, Double_t content) override
Definition GHSym.cxx:2842
TH1 * DrawCopy(Option_t *option="", const char *name_postfix="_copy") const override
Definition GHSym.cxx:2797
TH2F * GetMatrix(bool force=false)
Definition GHSym.cxx:2772
void Copy(TObject &rh) const override
Definition GHSym.cxx:2792
GHSymF()
Definition GHSym.cxx:2721
Definition GHSym.h:12
virtual Double_t DoIntegral(Int_t binx1, Int_t binx2, Int_t biny1, Int_t biny2, Double_t &error, Option_t *option, Bool_t doError=kFALSE) const
Definition GHSym.cxx:214
Double_t Integral(Option_t *option="") const override
Definition GHSym.cxx:1079
Long64_t Merge(TCollection *list) override
Definition GHSym.cxx:1448
virtual Double_t GetCovariance(Int_t axis1=1, Int_t axis2=2) const
Definition GHSym.cxx:922
Double_t fTsumwy
Total Sum of weight*Y.
Definition GHSym.h:114
void Copy(TObject &obj) const override
Definition GHSym.cxx:203
virtual Double_t GetBinWithContent2(Double_t c, Int_t &binx, Int_t &biny, Int_t firstxbin=1, Int_t lastxbin=-1, Int_t firstybin=1, Int_t lastybin=-1, Double_t maxdiff=0) const
Definition GHSym.cxx:834
virtual TProfile * Profile(const char *name="_pf", Int_t firstbin=1, Int_t lastbin=-1, Option_t *option="") const
Definition GHSym.cxx:1720
Double_t KolmogorovTest(const TH1 *h2, Option_t *option="") const override
Definition GHSym.cxx:1232
void SetCellError(Int_t binx, Int_t biny, Double_t content) override
Definition GHSym.cxx:2546
void FillRandom(const char *fname, Int_t ntimes=5000, TRandom *rng=nullptr) override
Definition GHSym.cxx:459
void GetStats(Double_t *stats) const override
Definition GHSym.cxx:996
void PutStats(Double_t *stats) override
Definition GHSym.cxx:2198
TH2 * fMatrix
! Transient pointer to the 2D-Matrix used in Draw() or GetMatrix()
Definition GHSym.h:117
Int_t BufferFill(Double_t, Double_t) override
Definition GHSym.h:26
void Reset(Option_t *option="") override
Definition GHSym.cxx:2516
void Smooth(Int_t ntimes=1, Option_t *option="") override
Definition GHSym.cxx:2602
TH1 * ShowBackground(Int_t niter=20, Option_t *option="same") override
Definition GHSym.cxx:2579
virtual Double_t IntegralAndError(Int_t firstxbin, Int_t lastxbin, Int_t firstybin, Int_t lastybin, Double_t &error, Option_t *option="") const
Definition GHSym.cxx:1100
TH2 * Matrix()
Definition GHSym.h:110
virtual Double_t GetCorrelationFactor(Int_t axis1=1, Int_t axis2=2) const
Definition GHSym.cxx:900
virtual void SetShowProjectionY(Int_t nbins=1)
Definition GHSym.cxx:2565
Double_t fTsumwy2
Total Sum of weight*Y*Y.
Definition GHSym.h:115
Int_t Fill(Double_t) override
Definition GHSym.cxx:281
Int_t BufferEmpty(Int_t action=0) override
Definition GHSym.cxx:77
Double_t Interpolate(Double_t) const override
Definition GHSym.cxx:1116
virtual TH1D * Projection(const char *name="_pr", Int_t firstBin=0, Int_t lastBin=-1, Option_t *option="") const
Definition GHSym.cxx:1957
Int_t ShowPeaks(Double_t sigma=2, Option_t *option="", Double_t threshold=0.05) override
Definition GHSym.cxx:2589
void FillN(Int_t, const Double_t *, const Double_t *, Int_t) override
Definition GHSym.h:34
Double_t GetCellError(Int_t binx, Int_t biny) const override
Definition GHSym.cxx:895
Int_t GetBin(Int_t binx, Int_t biny=0, Int_t binz=0) const override
Definition GHSym.cxx:813
virtual GHSym * Rebin2D(Int_t ngroup=2, const char *newname="")
Definition GHSym.cxx:2207
virtual void FitSlices(TF1 *f1=nullptr, Int_t firstbin=0, Int_t lastbin=-1, Int_t cut=0, Option_t *option="QNR", TObjArray *arr=nullptr)
Definition GHSym.cxx:640
Int_t FindFirstBinAbove(Double_t threshold=0, Int_t axis=1, Int_t firstBin=1, Int_t lastBin=-1) const override
Definition GHSym.cxx:574
Int_t FindLastBinAbove(Double_t threshold=0, Int_t axis=1, Int_t firstBin=1, Int_t lastBin=-1) const override
Definition GHSym.cxx:607
void SetCellContent(Int_t binx, Int_t biny, Double_t content) override
Definition GHSym.cxx:2541
Double_t GetCellContent(Int_t binx, Int_t biny) const override
Definition GHSym.cxx:890
virtual void GetRandom2(Double_t &x, Double_t &y)
Definition GHSym.cxx:953
virtual void SetShowProjectionX(Int_t nbins=1)
Definition GHSym.cxx:2551
GHSym()
Definition GHSym.cxx:30
Double_t fTsumwxy
Total Sum of weight*X*Y.
Definition GHSym.h:116