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