Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TProfile2Poly.cxx
Go to the documentation of this file.
1// @(#)root/hist:$Id$
2// Author: Filip Ilic
3
4/*************************************************************************
5 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12#include "TProfile2Poly.h"
13#include "TProfileHelper.h"
14#include "TMultiGraph.h"
15#include "TList.h"
16
17#include <iostream>
18#include <cassert>
19#include <cmath>
20#include <set>
21
22/** \class TProfile2Poly
23 \ingroup Histograms
242D Profile Histogram with Polygonal Bins.
25
26tprofile2polyRealisticModuleError.C and tprofile2polyRealistic.C illustrate how
27to use this class.
28*/
29
30/** \class TProfile2PolyBin
31 \ingroup Hist
32Helper class to represent a bin in the TProfile2Poly histogram
33*/
34
35
36////////////////////////////////////////////////////////////////////////////////
37/// TProfile2PolyBin constructor.
38
40{
41 fSumw = 0;
42 fSumvw = 0;
43 fSumw2 = 0;
44 fSumwv2 = 0;
45 fError = 0;
46 fAverage = 0;
48}
49
50////////////////////////////////////////////////////////////////////////////////
51/// TProfile2PolyBin constructor.
52
63
64////////////////////////////////////////////////////////////////////////////////
65/// Merge.
66
68{
69 this->fSumw += toMerge->fSumw;
70 this->fSumvw += toMerge->fSumvw;
71 this->fSumw2 += toMerge->fSumw2;
72 this->fSumwv2 += toMerge->fSumwv2;
73}
74
75////////////////////////////////////////////////////////////////////////////////
76/// Update.
77
79{
82 SetChanged(true);
83}
84
85////////////////////////////////////////////////////////////////////////////////
86/// Update average.
87
92
93////////////////////////////////////////////////////////////////////////////////
94/// Update error.
95
97{
98 Double_t tmp = 0;
99 if (fSumw != 0) tmp = std::sqrt((fSumwv2 / fSumw) - (fAverage * fAverage));
100
101 fError = tmp;
102
103 return;
104
105}
106
107////////////////////////////////////////////////////////////////////////////////
108/// Clear statistics.
109
111{
112 fSumw = 0;
113 fSumvw = 0;
114 fSumw2 = 0;
115 fSumwv2 = 0;
116 fError = 0;
117 fAverage = 0;
118}
119
120////////////////////////////////////////////////////////////////////////////////
121/// Fill.
122
124{
125 fSumw += weight;
126 fSumvw += value * weight;
127 fSumw2 += weight * weight;
128 fSumwv2 += weight * value * value;
129 this->Update();
130}
131
132////////////////////////////////////////////////////////////////////////////////
133/// TProfile2Poly constructor.
134
135TProfile2Poly::TProfile2Poly(const char *name, const char *title, Double_t xlow, Double_t xup, Double_t ylow,
136 Double_t yup)
137 : TH2Poly(name, title, xlow, xup, ylow, yup)
138{
139}
140
141////////////////////////////////////////////////////////////////////////////////
142/// TProfile2Poly constructor.
143
144TProfile2Poly::TProfile2Poly(const char *name, const char *title, Int_t nX, Double_t xlow, Double_t xup, Int_t nY,
145 Double_t ylow, Double_t yup)
146 : TH2Poly(name, title, nX, xlow, xup, nY, ylow, yup)
147{
148}
149
150////////////////////////////////////////////////////////////////////////////////
151/// Create bin.
152
154{
155 if (!poly) return nullptr;
156
157 if (fBins == nullptr) {
158 fBins = new TList();
159 fBins->SetOwner();
160 }
161
162 fNcells++;
164 return new TProfile2PolyBin(poly, ibin);
165}
166
167////////////////////////////////////////////////////////////////////////////////
168/// Fill
169
174
175////////////////////////////////////////////////////////////////////////////////
176/// Fill
177
179{
180 // Find region in which the hit occurred
182 if (tmp < 0) {
186 }
187
188 // Find the cell to which (x,y) coordinates belong to
189 Int_t n = (Int_t)(floor((xcoord - fXaxis.GetXmin()) / fStepX));
190 Int_t m = (Int_t)(floor((ycoord - fYaxis.GetXmin()) / fStepY));
191
192 // Make sure the array indices are correct.
193 if (n >= fCellX) n = fCellX - 1;
194 if (m >= fCellY) m = fCellY - 1;
195 if (n < 0) n = 0;
196 if (m < 0) m = 0;
197
198 // ------------ Update global (per histo) statistics
199 fTsumw += weight;
200 fTsumw2 += weight * weight;
201 fTsumwx += weight * xcoord;
202 fTsumwx2 += weight * xcoord * xcoord;
203 fTsumwy += weight * ycoord;
204 fTsumwy2 += weight * ycoord * ycoord;
205 fTsumwxy += weight * xcoord * ycoord;
206 fTsumwz += weight * value;
207 fTsumwz2 += weight * value * value;
208
209 // ------------ Update local (per bin) statistics
210 TProfile2PolyBin *bin;
211 TIter next(&fCells[n + fCellX * m]);
212 TObject *obj;
213 while ((obj = next())) {
214 bin = (TProfile2PolyBin *)obj;
215 if (bin->IsInside(xcoord, ycoord)) {
216 fEntries++;
217 bin->Fill(value, weight);
218 bin->Update();
219 bin->SetContent(bin->fAverage);
220 }
221 }
222
223 return tmp;
224}
225
226////////////////////////////////////////////////////////////////////////////////
227/// Merge
228
230{
231 Int_t size = in->GetSize();
232
233 std::vector<TProfile2Poly *> list;
234 list.reserve(size);
235
236 for (int i = 0; i < size; i++) {
237 list.push_back((TProfile2Poly *)((TList *)in)->At(i));
238 }
239 return this->Merge(list);
240}
241
242////////////////////////////////////////////////////////////////////////////////
243/// Merge
244
245Long64_t TProfile2Poly::Merge(const std::vector<TProfile2Poly *> &list)
246{
247 if (list.empty()) {
248 std::cout << "[FAIL] TProfile2Poly::Merge: No objects to be merged " << std::endl;
249 return -1;
250 }
251
252 // ------------ Check that bin numbers of TP2P's to be merged are equal
253 std::set<Int_t> numBinUnique;
254 for (const auto &histo : list) {
255 if (histo->fBins) numBinUnique.insert(histo->fBins->GetSize());
256 }
257 if (numBinUnique.size() != 1) {
258 std::cout << "[FAIL] TProfile2Poly::Merge: Bin numbers of TProfile2Polys to be merged differ!" << std::endl;
259 return -1;
260 }
261 Int_t nbins = *numBinUnique.begin();
262
263 // ------------ Update global (per histo) statistics
264 for (const auto &histo : list) {
265 this->fEntries += histo->fEntries;
266 this->fTsumw += histo->fTsumw;
267 this->fTsumw2 += histo->fTsumw2;
268 this->fTsumwx += histo->fTsumwx;
269 this->fTsumwx2 += histo->fTsumwx2;
270 this->fTsumwy += histo->fTsumwy;
271 this->fTsumwy2 += histo->fTsumwy2;
272 this->fTsumwxy += histo->fTsumwxy;
273 this->fTsumwz += histo->fTsumwz;
274 this->fTsumwz2 += histo->fTsumwz2;
275
276 // Merge overflow bins
277 for (Int_t i = 0; i < kNOverflow; ++i) {
278 this->fOverflowBins[i].Merge(&histo->fOverflowBins[i]);
279 }
280 }
281
282 // ------------ Update local (per bin) statistics
283 TProfile2PolyBin *dst = nullptr;
284 TProfile2PolyBin *src = nullptr;
285 for (Int_t i = 0; i < nbins; i++) {
286 dst = (TProfile2PolyBin *)fBins->At(i);
287
288 for (const auto &e : list) {
289 src = (TProfile2PolyBin *)e->fBins->At(i);
290 dst->Merge(src);
291 }
292
293 dst->Update();
294 }
295
296 this->SetContentToAverage();
297 return 1;
298}
299
300////////////////////////////////////////////////////////////////////////////////
301/// Set content to average.
302
304{
305 Int_t nbins = fBins ? fBins->GetSize() : 0;
306 for (Int_t i = 0; i < nbins; i++) {
308 bin->Update();
309 bin->SetContent(bin->fAverage);
310 }
311 for (Int_t i = 0; i < kNOverflow; ++i) {
313 bin.Update();
314 bin.SetContent(bin.fAverage);
315 }
316}
317
318////////////////////////////////////////////////////////////////////////////////
319/// Set content to error.
320
322{
323 Int_t nbins = fBins ? fBins->GetSize() : 0;
324 for (Int_t i = 0; i < nbins; i++) {
326 bin->Update();
327 bin->SetContent(bin->fError);
328 }
329 for (Int_t i = 0; i < kNOverflow; ++i) {
331 bin.Update();
332 bin.SetContent(bin.fError);
333 }
334}
335
336////////////////////////////////////////////////////////////////////////////////
337/// Get bin content.
338
340{
341 if (bin > GetNumberOfBins() || bin == 0 || bin < -kNOverflow) return 0;
342 if (bin<0) return fOverflowBins[-bin - 1].GetContent();
343 return ((TProfile2PolyBin*) fBins->At(bin-1))->GetContent();
344}
345
346
347////////////////////////////////////////////////////////////////////////////////
348/// Get bin effective entries.
349
351{
352 if (bin > GetNumberOfBins() || bin == 0 || bin < -kNOverflow) return 0;
353 if (bin < 0) return fOverflowBins[-bin - 1].GetEffectiveEntries();
354 return ((TProfile2PolyBin *)fBins->At(bin - 1))->GetEffectiveEntries();
355}
356
357////////////////////////////////////////////////////////////////////////////////
358/// Get bin entries.
359
361{
362 if (bin > GetNumberOfBins() || bin == 0 || bin < -kNOverflow) return 0;
363 if (bin < 0) return fOverflowBins[-bin - 1].GetEntries();
364 return ((TProfile2PolyBin *)fBins->At(bin - 1))->GetEntries();
365}
366
367////////////////////////////////////////////////////////////////////////////////
368/// Get bin entries W2.
369
371{
372 if (bin > GetNumberOfBins() || bin == 0 || bin < -kNOverflow) return 0;
373 if (bin < 0) return fOverflowBins[-bin - 1].GetEntriesW2();
374 return ((TProfile2PolyBin *)fBins->At(bin - 1))->GetEntriesW2();
375}
376
377////////////////////////////////////////////////////////////////////////////////
378/// Get bin entries VW.
379
381{
382 if (bin > GetNumberOfBins() || bin == 0 || bin < -kNOverflow) return 0;
383 if (bin < 0) return fOverflowBins[-bin - 1].GetEntriesVW();
384 return ((TProfile2PolyBin *)fBins->At(bin - 1))->GetEntriesVW();
385}
386
387////////////////////////////////////////////////////////////////////////////////
388/// Get bin entries WV2.
389
391{
392 if (bin > GetNumberOfBins() || bin == 0 || bin < -kNOverflow) return 0;
393 if (bin < 0) return fOverflowBins[-bin - 1].GetEntriesWV2();
394 return ((TProfile2PolyBin *)fBins->At(bin - 1))->GetEntriesWV2();
395}
396
397////////////////////////////////////////////////////////////////////////////////
398/// Get bin error.
399
401{
402 Double_t tmp = 0;
403 if (bin > GetNumberOfBins() || bin == 0 || bin < -kNOverflow) return 0;
404 if (bin < 0)
405 tmp = fOverflowBins[-bin - 1].GetError();
406 else
407 tmp = ((TProfile2PolyBin *)fBins->At(bin - 1))->GetError();
408
409 return (fErrorMode == kERRORSPREAD) ? tmp : tmp / std::sqrt(GetBinEffectiveEntries(bin));
410
411}
412
413////////////////////////////////////////////////////////////////////////////////
414/// Fill the array stats from the contents of this profile.
415/// The array stats must be correctly dimensioned in the calling program.
416///
417/// - stats[0] = sumw
418/// - stats[1] = sumw2
419/// - stats[2] = sumwx
420/// - stats[3] = sumwx2
421/// - stats[4] = sumwy
422/// - stats[5] = sumwy2
423/// - stats[6] = sumwxy
424/// - stats[7] = sumwz
425/// - stats[8] = sumwz2
426///
427/// If no axis-subrange is specified (via TAxis::SetRange), the array stats
428/// is simply a copy of the statistics quantities computed at filling time.
429/// If a sub-range is specified, the function recomputes these quantities
430/// from the bin contents in the current axis range.
431
433{
434 stats[0] = fTsumw;
435 stats[1] = fTsumw2;
436 stats[2] = fTsumwx;
437 stats[3] = fTsumwx2;
438 stats[4] = fTsumwy;
439 stats[5] = fTsumwy2;
440 stats[6] = fTsumwxy;
441 stats[7] = fTsumwz;
442 stats[8] = fTsumwz2;
443}
444
445////////////////////////////////////////////////////////////////////////////////
446/// Print overflow regions.
447
449{
450 Double_t total = 0;
451 Double_t cont = 0;
452 for (Int_t i = 0; i < kNOverflow; ++i) {
454 total += cont;
455 std::cout << "\t" << cont << "\t";
456 if ((i + 1) % 3 == 0) std::cout << std::endl;
457 }
458
459 std::cout << "Total: " << total << std::endl;
460}
461
462////////////////////////////////////////////////////////////////////////////////
463/// Reset
464
466{
467 TIter next(fBins);
468 TObject *obj;
469 TProfile2PolyBin *bin;
470
471 // Clears bin contents
472 while ((obj = next())) {
473 bin = (TProfile2PolyBin *)obj;
474 bin->ClearContent();
475 bin->ClearStats();
476 }
477 TH2::Reset(opt);
478}
479
480////////////////////////////////////////////////////////////////////////////////
481/// The overflow regions are calculated by considering x, y coordinates.
482/// The Middle bin at -5 contains all the TProfile2Poly bins.
483///
484/// ~~~ {.cpp}
485/// -0 -1 -2
486/// ________
487/// -1: |__|__|__|
488/// -4: |__|__|__|
489/// -7: |__|__|__|
490/// ~~~
491
493{
494
495
496 Int_t region = 0;
497
498 if (fNcells <= kNOverflow) return 0;
499
500 // --- y offset
501 if (y > fYaxis.GetXmax())
502 region += -1;
503 else if (y > fYaxis.GetXmin())
504 region += -4;
505 else
506 region += -7;
507
508 // --- x offset
509 if (x > fXaxis.GetXmax())
510 region += -2;
511 else if (x > fXaxis.GetXmin())
512 region += -1;
513 else
514 region += 0;
515
516 return region;
517}
518
519////////////////////////////////////////////////////////////////////////////////
520/// Set error option.
521
#define e(i)
Definition RSha256.hxx:103
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
static unsigned int total
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t src
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
char name[80]
Definition TGX11.cxx:110
EErrorType
Definition TProfile.h:28
@ kERRORSPREAD
Definition TProfile.h:28
@ kERRORMEAN
Definition TProfile.h:28
const_iterator begin() const
Double_t GetXmax() const
Definition TAxis.h:142
Double_t GetXmin() const
Definition TAxis.h:141
Collection abstract base class.
Definition TCollection.h:65
virtual void SetOwner(Bool_t enable=kTRUE)
Set whether this collection is the owner (enable==true) of its content.
virtual Int_t GetSize() const
Return the capacity of the collection, i.e.
Int_t fNcells
Number of bins(1D), cells (2D) +U/Overflows.
Definition TH1.h:150
Double_t fTsumw
Total Sum of weights.
Definition TH1.h:157
Double_t fTsumw2
Total Sum of squares of weights.
Definition TH1.h:158
Double_t fTsumwx2
Total Sum of weight*X*X.
Definition TH1.h:160
Double_t fEntries
Number of entries.
Definition TH1.h:156
TAxis fXaxis
X axis descriptor.
Definition TH1.h:151
TAxis fYaxis
Y axis descriptor.
Definition TH1.h:152
Double_t fTsumwx
Total Sum of weight*X.
Definition TH1.h:159
Helper class to represent a bin in the TH2Poly histogram.
Definition TH2Poly.h:25
void ClearContent()
Definition TH2Poly.h:32
Bool_t IsInside(Double_t x, Double_t y) const
Return "true" if the point (x,y) is inside the bin.
Definition TH2Poly.cxx:1643
void SetContent(Double_t content)
Definition TH2Poly.h:45
Double_t GetContent() const
Definition TH2Poly.h:35
void SetChanged(Bool_t flag)
Definition TH2Poly.h:44
2D Histogram with Polygonal Bins
Definition TH2Poly.h:66
@ kNOverflow
Definition TH2Poly.h:159
TList * fBins
List of bins. The list owns the contained objects.
Definition TH2Poly.h:172
Int_t GetNumberOfBins() const
Return the number of bins : it should be the size of the bin list.
Definition TH2Poly.cxx:863
Double_t fStepX
Definition TH2Poly.h:166
Int_t fCellX
Number of partition cells in the x-direction of the histogram.
Definition TH2Poly.h:162
Double_t fStepY
Dimensions of a partition cell.
Definition TH2Poly.h:166
TList * fCells
[fNCells] The array of TLists that store the bins that intersect with each cell. List do not own the ...
Definition TH2Poly.h:165
Int_t fCellY
Number of partition cells in the y-direction of the histogram.
Definition TH2Poly.h:163
Double_t fTsumwxy
Total Sum of weight*X*Y.
Definition TH2.h:36
void Reset(Option_t *option="") override
Reset this histogram: contents, errors, etc.
Definition TH2.cxx:2562
Double_t fTsumwy2
Total Sum of weight*Y*Y.
Definition TH2.h:35
Double_t fTsumwy
Total Sum of weight*Y.
Definition TH2.h:34
A doubly linked list.
Definition TList.h:38
TObject * At(Int_t idx) const override
Returns the object at position idx. Returns 0 if idx is out of range.
Definition TList.cxx:354
Mother of all ROOT objects.
Definition TObject.h:41
Helper class to represent a bin in the TProfile2Poly histogram.
Double_t GetEntriesW2() const
void Update()
Update.
Double_t GetEntries() const
void Fill(Double_t value, Double_t weight)
Fill.
Double_t GetError() const
Double_t GetEffectiveEntries() const
TProfile2PolyBin()
TProfile2PolyBin constructor.
void UpdateAverage()
Update average.
Double_t GetEntriesVW() const
void UpdateError()
Update error.
EErrorType fErrorMode
Double_t GetEntriesWV2() const
void Merge(const TProfile2PolyBin *toMerge)
Merge.
void ClearStats()
Clear statistics.
2D Profile Histogram with Polygonal Bins.
void PrintOverflowRegions()
Print overflow regions.
Double_t GetOverflowContent(Int_t idx)
EErrorType fErrorMode
Int_t Fill(Double_t xcoord, Double_t ycoord, Double_t value) override
Fill.
Long64_t Merge(const std::vector< TProfile2Poly * > &list)
Merge.
Double_t fTsumwz2
Double_t GetBinEffectiveEntries(Int_t bin) const
Get bin effective entries.
Double_t GetBinError(Int_t bin) const override
Get bin error.
Double_t GetBinEntriesVW(Int_t bin) const
Get bin entries VW.
Double_t GetBinEntriesW2(Int_t bin) const
Get bin entries W2.
friend class TProfile2PolyBin
void SetContentToError()
Set content to error.
Double_t GetBinContent(Int_t bin) const override
Get bin content.
Int_t GetOverflowRegionFromCoordinates(Double_t x, Double_t y)
The overflow regions are calculated by considering x, y coordinates.
void GetStats(Double_t *stats) const override
Fill the array stats from the contents of this profile.
void Reset(Option_t *option="") override
Reset.
void SetContentToAverage()
Set content to average.
Double_t GetBinEntriesWV2(Int_t bin) const
Get bin entries WV2.
Double_t GetBinEntries(Int_t bin) const
Get bin entries.
void SetErrorOption(EErrorType type)
Set error option.
Int_t OverflowIdxToArrayIdx(Int_t val)
TProfile2PolyBin fOverflowBins[kNOverflow]
TProfile2PolyBin * CreateBin(TObject *poly) override
Create bin.
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
TMarker m
Definition textangle.C:8