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
36
37////////////////////////////////////////////////////////////////////////////////
38/// TProfile2PolyBin constructor.
39
41{
42 fSumw = 0;
43 fSumvw = 0;
44 fSumw2 = 0;
45 fSumwv2 = 0;
46 fError = 0;
47 fAverage = 0;
49}
50
51////////////////////////////////////////////////////////////////////////////////
52/// TProfile2PolyBin constructor.
53
54TProfile2PolyBin::TProfile2PolyBin(TObject *poly, Int_t bin_number) : TH2PolyBin(poly, bin_number)
55{
56 fSumw = 0;
57 fSumvw = 0;
58 fSumw2 = 0;
59 fSumwv2 = 0;
60 fError = 0;
61 fAverage = 0;
63}
64
65////////////////////////////////////////////////////////////////////////////////
66/// Merge.
67
69{
70 this->fSumw += toMerge->fSumw;
71 this->fSumvw += toMerge->fSumvw;
72 this->fSumw2 += toMerge->fSumw2;
73 this->fSumwv2 += toMerge->fSumwv2;
74}
75
76////////////////////////////////////////////////////////////////////////////////
77/// Update.
78
80{
83 SetChanged(true);
84}
85
86////////////////////////////////////////////////////////////////////////////////
87/// Update average.
88
90{
91 if (fSumw != 0) fAverage = fSumvw / fSumw;
92}
93
94////////////////////////////////////////////////////////////////////////////////
95/// Update error.
96
98{
99 Double_t tmp = 0;
100 if (fSumw != 0) tmp = std::sqrt((fSumwv2 / fSumw) - (fAverage * fAverage));
101
102 fError = tmp;
103
104 return;
105
106}
107
108////////////////////////////////////////////////////////////////////////////////
109/// Clear statistics.
110
112{
113 fSumw = 0;
114 fSumvw = 0;
115 fSumw2 = 0;
116 fSumwv2 = 0;
117 fError = 0;
118 fAverage = 0;
119}
120
121////////////////////////////////////////////////////////////////////////////////
122/// Fill.
123
125{
126 fSumw += weight;
127 fSumvw += value * weight;
128 fSumw2 += weight * weight;
129 fSumwv2 += weight * value * value;
130 this->Update();
131}
132
133////////////////////////////////////////////////////////////////////////////////
134/// TProfile2Poly constructor.
135
136TProfile2Poly::TProfile2Poly(const char *name, const char *title, Double_t xlow, Double_t xup, Double_t ylow,
137 Double_t yup)
138 : TH2Poly(name, title, xlow, xup, ylow, yup)
139{
140}
141
142////////////////////////////////////////////////////////////////////////////////
143/// TProfile2Poly constructor.
144
145TProfile2Poly::TProfile2Poly(const char *name, const char *title, Int_t nX, Double_t xlow, Double_t xup, Int_t nY,
146 Double_t ylow, Double_t yup)
147 : TH2Poly(name, title, nX, xlow, xup, nY, ylow, yup)
148{
149}
150
151////////////////////////////////////////////////////////////////////////////////
152/// Create bin.
153
155{
156 if (!poly) return nullptr;
157
158 if (fBins == nullptr) {
159 fBins = new TList();
160 fBins->SetOwner();
161 }
162
163 fNcells++;
164 Int_t ibin = fNcells - kNOverflow;
165 return new TProfile2PolyBin(poly, ibin);
166}
167
168////////////////////////////////////////////////////////////////////////////////
169/// Fill
170
172{
173 return Fill(xcoord, ycoord, value, 1);
174}
175
176////////////////////////////////////////////////////////////////////////////////
177/// Fill
178
180{
181 // Find region in which the hit occurred
182 Int_t tmp = GetOverflowRegionFromCoordinates(xcoord, ycoord);
183 if (tmp < 0) {
184 Int_t overflow_idx = OverflowIdxToArrayIdx(tmp);
185 fOverflowBins[overflow_idx].Fill(value, weight);
186 fOverflowBins[overflow_idx].SetContent(fOverflowBins[overflow_idx].fAverage );
187 }
188
189 // Find the cell to which (x,y) coordinates belong to
190 Int_t n = (Int_t)(floor((xcoord - fXaxis.GetXmin()) / fStepX));
191 Int_t m = (Int_t)(floor((ycoord - fYaxis.GetXmin()) / fStepY));
192
193 // Make sure the array indices are correct.
194 if (n >= fCellX) n = fCellX - 1;
195 if (m >= fCellY) m = fCellY - 1;
196 if (n < 0) n = 0;
197 if (m < 0) m = 0;
198
199 // ------------ Update global (per histo) statistics
200 fTsumw += weight;
201 fTsumw2 += weight * weight;
202 fTsumwx += weight * xcoord;
203 fTsumwx2 += weight * xcoord * xcoord;
204 fTsumwy += weight * ycoord;
205 fTsumwy2 += weight * ycoord * ycoord;
206 fTsumwxy += weight * xcoord * ycoord;
207 fTsumwz += weight * value;
208 fTsumwz2 += weight * value * value;
209
210 // ------------ Update local (per bin) statistics
211 TProfile2PolyBin *bin;
212 TIter next(&fCells[n + fCellX * m]);
213 TObject *obj;
214 while ((obj = next())) {
215 bin = (TProfile2PolyBin *)obj;
216 if (bin->IsInside(xcoord, ycoord)) {
217 fEntries++;
218 bin->Fill(value, weight);
219 bin->Update();
220 bin->SetContent(bin->fAverage);
221 }
222 }
223
224 return tmp;
225}
226
227////////////////////////////////////////////////////////////////////////////////
228/// Merge
229
231{
232 Int_t size = in->GetSize();
233
234 std::vector<TProfile2Poly *> list;
235 list.reserve(size);
236
237 for (int i = 0; i < size; i++) {
238 list.push_back((TProfile2Poly *)((TList *)in)->At(i));
239 }
240 return this->Merge(list);
241}
242
243////////////////////////////////////////////////////////////////////////////////
244/// Merge
245
246Long64_t TProfile2Poly::Merge(const std::vector<TProfile2Poly *> &list)
247{
248 if (list.empty()) {
249 std::cout << "[FAIL] TProfile2Poly::Merge: No objects to be merged " << std::endl;
250 return -1;
251 }
252
253 // ------------ Check that bin numbers of TP2P's to be merged are equal
254 std::set<Int_t> numBinUnique;
255 for (const auto &histo : list) {
256 if (histo->fBins) numBinUnique.insert(histo->fBins->GetSize());
257 }
258 if (numBinUnique.size() != 1) {
259 std::cout << "[FAIL] TProfile2Poly::Merge: Bin numbers of TProfile2Polys to be merged differ!" << std::endl;
260 return -1;
261 }
262 Int_t nbins = *numBinUnique.begin();
263
264 // ------------ Update global (per histo) statistics
265 for (const auto &histo : list) {
266 this->fEntries += histo->fEntries;
267 this->fTsumw += histo->fTsumw;
268 this->fTsumw2 += histo->fTsumw2;
269 this->fTsumwx += histo->fTsumwx;
270 this->fTsumwx2 += histo->fTsumwx2;
271 this->fTsumwy += histo->fTsumwy;
272 this->fTsumwy2 += histo->fTsumwy2;
273 this->fTsumwxy += histo->fTsumwxy;
274 this->fTsumwz += histo->fTsumwz;
275 this->fTsumwz2 += histo->fTsumwz2;
276
277 // Merge overflow bins
278 for (Int_t i = 0; i < kNOverflow; ++i) {
279 this->fOverflowBins[i].Merge(&histo->fOverflowBins[i]);
280 }
281 }
282
283 // ------------ Update local (per bin) statistics
284 TProfile2PolyBin *dst = nullptr;
285 TProfile2PolyBin *src = nullptr;
286 for (Int_t i = 0; i < nbins; i++) {
287 dst = (TProfile2PolyBin *)fBins->At(i);
288
289 for (const auto &e : list) {
290 src = (TProfile2PolyBin *)e->fBins->At(i);
291 dst->Merge(src);
292 }
293
294 dst->Update();
295 }
296
297 this->SetContentToAverage();
298 return 1;
299}
300
301////////////////////////////////////////////////////////////////////////////////
302/// Set content to average.
303
305{
306 Int_t nbins = fBins ? fBins->GetSize() : 0;
307 for (Int_t i = 0; i < nbins; i++) {
309 bin->Update();
310 bin->SetContent(bin->fAverage);
311 }
312 for (Int_t i = 0; i < kNOverflow; ++i) {
314 bin.Update();
315 bin.SetContent(bin.fAverage);
316 }
317}
318
319////////////////////////////////////////////////////////////////////////////////
320/// Set content to error.
321
323{
324 Int_t nbins = fBins ? fBins->GetSize() : 0;
325 for (Int_t i = 0; i < nbins; i++) {
327 bin->Update();
328 bin->SetContent(bin->fError);
329 }
330 for (Int_t i = 0; i < kNOverflow; ++i) {
332 bin.Update();
333 bin.SetContent(bin.fError);
334 }
335}
336
337////////////////////////////////////////////////////////////////////////////////
338/// Get bin content.
339
341{
342 if (bin > GetNumberOfBins() || bin == 0 || bin < -kNOverflow) return 0;
343 if (bin<0) return fOverflowBins[-bin - 1].GetContent();
344 return ((TProfile2PolyBin*) fBins->At(bin-1))->GetContent();
345}
346
347
348////////////////////////////////////////////////////////////////////////////////
349/// Get bin effective entries.
350
352{
353 if (bin > GetNumberOfBins() || bin == 0 || bin < -kNOverflow) return 0;
354 if (bin < 0) return fOverflowBins[-bin - 1].GetEffectiveEntries();
355 return ((TProfile2PolyBin *)fBins->At(bin - 1))->GetEffectiveEntries();
356}
357
358////////////////////////////////////////////////////////////////////////////////
359/// Get bin entries.
360
362{
363 if (bin > GetNumberOfBins() || bin == 0 || bin < -kNOverflow) return 0;
364 if (bin < 0) return fOverflowBins[-bin - 1].GetEntries();
365 return ((TProfile2PolyBin *)fBins->At(bin - 1))->GetEntries();
366}
367
368////////////////////////////////////////////////////////////////////////////////
369/// Get bin entries W2.
370
372{
373 if (bin > GetNumberOfBins() || bin == 0 || bin < -kNOverflow) return 0;
374 if (bin < 0) return fOverflowBins[-bin - 1].GetEntriesW2();
375 return ((TProfile2PolyBin *)fBins->At(bin - 1))->GetEntriesW2();
376}
377
378////////////////////////////////////////////////////////////////////////////////
379/// Get bin entries VW.
380
382{
383 if (bin > GetNumberOfBins() || bin == 0 || bin < -kNOverflow) return 0;
384 if (bin < 0) return fOverflowBins[-bin - 1].GetEntriesVW();
385 return ((TProfile2PolyBin *)fBins->At(bin - 1))->GetEntriesVW();
386}
387
388////////////////////////////////////////////////////////////////////////////////
389/// Get bin entries WV2.
390
392{
393 if (bin > GetNumberOfBins() || bin == 0 || bin < -kNOverflow) return 0;
394 if (bin < 0) return fOverflowBins[-bin - 1].GetEntriesWV2();
395 return ((TProfile2PolyBin *)fBins->At(bin - 1))->GetEntriesWV2();
396}
397
398////////////////////////////////////////////////////////////////////////////////
399/// Get bin error.
400
402{
403 Double_t tmp = 0;
404 if (bin > GetNumberOfBins() || bin == 0 || bin < -kNOverflow) return 0;
405 if (bin < 0)
406 tmp = fOverflowBins[-bin - 1].GetError();
407 else
408 tmp = ((TProfile2PolyBin *)fBins->At(bin - 1))->GetError();
409
410 return (fErrorMode == kERRORSPREAD) ? tmp : tmp / std::sqrt(GetBinEffectiveEntries(bin));
411
412}
413
414////////////////////////////////////////////////////////////////////////////////
415/// Fill the array stats from the contents of this profile.
416/// The array stats must be correctly dimensioned in the calling program.
417///
418/// - stats[0] = sumw
419/// - stats[1] = sumw2
420/// - stats[2] = sumwx
421/// - stats[3] = sumwx2
422/// - stats[4] = sumwy
423/// - stats[5] = sumwy2
424/// - stats[6] = sumwxy
425/// - stats[7] = sumwz
426/// - stats[8] = sumwz2
427///
428/// If no axis-subrange is specified (via TAxis::SetRange), the array stats
429/// is simply a copy of the statistics quantities computed at filling time.
430/// If a sub-range is specified, the function recomputes these quantities
431/// from the bin contents in the current axis range.
432
434{
435 stats[0] = fTsumw;
436 stats[1] = fTsumw2;
437 stats[2] = fTsumwx;
438 stats[3] = fTsumwx2;
439 stats[4] = fTsumwy;
440 stats[5] = fTsumwy2;
441 stats[6] = fTsumwxy;
442 stats[7] = fTsumwz;
443 stats[8] = fTsumwz2;
444}
445
446////////////////////////////////////////////////////////////////////////////////
447/// Print overflow regions.
448
450{
451 Double_t total = 0;
452 Double_t cont = 0;
453 for (Int_t i = 0; i < kNOverflow; ++i) {
454 cont = GetOverflowContent(i);
455 total += cont;
456 std::cout << "\t" << cont << "\t";
457 if ((i + 1) % 3 == 0) std::cout << std::endl;
458 }
459
460 std::cout << "Total: " << total << std::endl;
461}
462
463////////////////////////////////////////////////////////////////////////////////
464/// Reset
465
467{
468 TIter next(fBins);
469 TObject *obj;
470 TProfile2PolyBin *bin;
471
472 // Clears bin contents
473 while ((obj = next())) {
474 bin = (TProfile2PolyBin *)obj;
475 bin->ClearContent();
476 bin->ClearStats();
477 }
478 TH2::Reset(opt);
479}
480
481////////////////////////////////////////////////////////////////////////////////
482/// The overflow regions are calculated by considering x, y coordinates.
483/// The Middle bin at -5 contains all the TProfile2Poly bins.
484///
485/// ~~~ {.cpp}
486/// -0 -1 -2
487/// ________
488/// -1: |__|__|__|
489/// -4: |__|__|__|
490/// -7: |__|__|__|
491/// ~~~
492
494{
495
496
497 Int_t region = 0;
498
499 if (fNcells <= kNOverflow) return 0;
500
501 // --- y offset
502 if (y > fYaxis.GetXmax())
503 region += -1;
504 else if (y > fYaxis.GetXmin())
505 region += -4;
506 else
507 region += -7;
508
509 // --- x offset
510 if (x > fXaxis.GetXmax())
511 region += -2;
512 else if (x > fXaxis.GetXmin())
513 region += -1;
514 else
515 region += 0;
516
517 return region;
518}
519
520////////////////////////////////////////////////////////////////////////////////
521/// Set error option.
522
524{
526}
#define e(i)
Definition RSha256.hxx:103
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
int Int_t
Definition RtypesCore.h:45
long long Long64_t
Definition RtypesCore.h:80
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:377
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
Double_t GetXmax() const
Definition TAxis.h:140
Double_t GetXmin() const
Definition TAxis.h:139
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:89
Double_t fTsumw
Total Sum of weights.
Definition TH1.h:96
Double_t fTsumw2
Total Sum of squares of weights.
Definition TH1.h:97
Double_t fTsumwx2
Total Sum of weight*X*X.
Definition TH1.h:99
Double_t fEntries
Number of entries.
Definition TH1.h:95
TAxis fXaxis
X axis descriptor.
Definition TH1.h:90
TAxis fYaxis
Y axis descriptor.
Definition TH1.h:91
Double_t fTsumwx
Total Sum of weight*X.
Definition TH1.h:98
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:1661
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:157
TList * fBins
List of bins. The list owns the contained objects.
Definition TH2Poly.h:170
Int_t GetNumberOfBins() const
Return the number of bins : it should be the size of the bin list.
Definition TH2Poly.cxx:865
Double_t fStepX
Definition TH2Poly.h:164
Int_t fCellX
Number of partition cells in the x-direction of the histogram.
Definition TH2Poly.h:160
Double_t fStepY
Dimensions of a partition cell.
Definition TH2Poly.h:164
TList * fCells
[fNCells] The array of TLists that store the bins that intersect with each cell. List do not own the ...
Definition TH2Poly.h:163
Int_t fCellY
Number of partition cells in the y-direction of the histogram.
Definition TH2Poly.h:161
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:2600
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:355
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