Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
TH2Poly.cxx
Go to the documentation of this file.
1 // @(#)root/hist:$Id$
2 // TH2Poly v2.1
3 // Author: Olivier Couet, Deniz Gunceler, Danilo Piparo
4 
5 /*************************************************************************
6  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
7  * All rights reserved. *
8  * *
9  * For the licensing terms see $ROOTSYS/LICENSE. *
10  * For the list of contributors see $ROOTSYS/README/CREDITS. *
11  *************************************************************************/
12 
13 #include "TH2Poly.h"
14 #include "TMultiGraph.h"
15 #include "TGraph.h"
16 #include "TClass.h"
17 #include "TList.h"
18 #include "TMath.h"
19 #include <cassert>
20 
21 ClassImp(TH2Poly);
22 
23 /** \class TH2Poly
24  \ingroup Hist
25 2D Histogram with Polygonal Bins
26 
27 ## Overview
28 `TH2Poly` is a 2D Histogram class (TH2) allowing to define polygonal
29 bins of arbitrary shape.
30 
31 Each bin in the `TH2Poly` histogram is a `TH2PolyBin` object.
32 `TH2PolyBin` is a very simple class containing the vertices (stored
33 as `TGraph`s or `TMultiGraph`s ) and contents of the polygonal
34 bin as well as several related functions.
35 
36 Essentially, a `TH2Poly` is a TList of `TH2PolyBin` objects
37 with methods to manipulate them.
38 
39 Bins are defined using one of the `AddBin()` methods. The bin definition
40 should be done before filling.
41 
42 The histogram can be filled with `Fill(Double_t x, Double_t y, Double_t w)
43 `. `w` is the weight.
44 If no weight is specified, it is assumed to be 1.
45 
46 Not all histogram's area need to be binned. Filling an area without bins,
47 will falls into the overflows. Adding a bin is not retroactive; it doesn't
48 affect previous fillings. A `Fill()` call, that
49 was previously ignored due to the lack of a bin at the specified location, is
50 not reconsidered when that location is binned later.
51 
52 If there are two overlapping bins, the first one in the list will be incremented
53 by `Fill()`.
54 
55 The histogram may automatically extends its limits if a bin outside the
56 histogram limits is added. This is done when the default constructor (with no
57 arguments) is used. It generates a histogram with no limits along the X and Y
58 axis. Adding bins to it will extend it up to a proper size.
59 
60 `TH2Poly` implements a partitioning algorithm to speed up bins' filling
61 (see the "Partitioning Algorithm" section for details).
62 The partitioning algorithm divides the histogram into regions called cells.
63 The bins that each cell intersects are recorded in an array of `TList`s.
64 When a coordinate in the histogram is to be filled; the method (quickly) finds
65 which cell the coordinate belongs. It then only loops over the bins
66 intersecting that cell to find the bin the input coordinate corresponds to.
67 The partitioning of the histogram is updated continuously as each bin is added.
68 The default number of cells on each axis is 25. This number could be set to
69 another value in the constructor or adjusted later by calling the
70 `ChangePartition(Int_t, Int_t)` method. The partitioning algorithm is
71 considerably faster than the brute force algorithm (i.e. checking if each bin
72 contains the input coordinates), especially if the histogram is to be filled
73 many times.
74 
75 The following very simple macro shows how to build and fill a `TH2Poly`:
76 ~~~ {.cpp}
77 {
78  TH2Poly *h2p = new TH2Poly();
79 
80  Double_t x1[] = {0, 5, 6};
81  Double_t y1[] = {0, 0, 5};
82  Double_t x2[] = {0, -1, -1, 0};
83  Double_t y2[] = {0, 0, -1, 3};
84  Double_t x3[] = {4, 3, 0, 1, 2.4};
85  Double_t y3[] = {4, 3.7, 1, 3.7, 2.5};
86 
87  h2p->AddBin(3, x1, y1);
88  h2p->AddBin(4, x2, y2);
89  h2p->AddBin(5, x3, y3);
90 
91  h2p->Fill(0.1, 0.01, 3);
92  h2p->Fill(-0.5, -0.5, 7);
93  h2p->Fill(-0.7, -0.5, 1);
94  h2p->Fill(1, 3, 1.5);
95 }
96 ~~~
97 
98 More examples can be found in th2polyBoxes.C, th2polyEurope.C, th2polyHoneycomb.C
99 and th2polyUSA.C.
100 
101 ## Partitioning Algorithm
102 The partitioning algorithm forms an essential part of the `TH2Poly`
103 class. It is implemented to speed up the filling of bins.
104 
105 With the brute force approach, the filling is done in the following way: An
106 iterator loops over all bins in the `TH2Poly` and invokes the
107 method `IsInside()` for each of them.
108 This method checks if the input location is in that bin. If the filling
109 coordinate is inside, the bin is filled. Looping over all the bin is
110 very slow.
111 
112 The alternative is to divide the histogram into virtual rectangular regions
113 called "cells". Each cell stores the pointers of the bins intersecting it.
114 When a coordinate is to be filled, the method finds which cell the coordinate
115 falls into. Since the cells are rectangular, this can be done very quickly.
116 It then only loops over the bins associated with that cell and calls `IsInside()`
117 only on that bins. This reduces considerably the number of bins on which `IsInside()`
118 is called and therefore speed up by a huge factor the filling compare to the brute force
119 approach where `IsInside()` is called for all bins.
120 
121 The addition of bins to the appropriate cells is done when the bin is added
122 to the histogram. To do this, `AddBin()` calls the
123 `AddBinToPartition()` method.
124 This method adds the input bin to the partitioning matrix.
125 
126 The number of partition cells per axis can be specified in the constructor.
127 If it is not specified, the default value of 25 along each axis will be
128 assigned. This value was chosen because it is small enough to avoid slowing
129 down AddBin(), while being large enough to enhance Fill() by a considerable
130 amount. Regardless of how it is initialized at construction time, it can be
131 changed later with the `ChangePartition()` method.
132 `ChangePartition()` deletes the
133 old partition matrix and generates a new one with the specified number of cells
134 on each axis.
135 
136 The optimum number of partition cells per axis changes with the number of
137 times `Fill()` will be called. Although partitioning greatly speeds up
138 filling, it also adds a constant time delay into the code. When `Fill()`
139 is to be called many times, it is more efficient to divide the histogram into
140 a large number cells. However, if the histogram is to be filled only a few
141 times, it is better to divide into a small number of cells.
142 */
143 
144 ////////////////////////////////////////////////////////////////////////////////
145 /// Default Constructor. No boundaries specified.
146 
147 TH2Poly::TH2Poly()
148 {
149  Initialize(0., 0., 0., 0., 25, 25);
150  SetName("NoName");
151  SetTitle("NoTitle");
152  SetFloat();
153 }
154 
155 ////////////////////////////////////////////////////////////////////////////////
156 /// Constructor with specified name and boundaries,
157 /// but no partition cell number.
158 
159 TH2Poly::TH2Poly(const char *name,const char *title, Double_t xlow,Double_t xup
160  , Double_t ylow,Double_t yup)
161 {
162  Initialize(xlow, xup, ylow, yup, 25, 25);
163  SetName(name);
164  SetTitle(title);
165  SetFloat(kFALSE);
166 }
167 
168 ////////////////////////////////////////////////////////////////////////////////
169 /// Constructor with specified name and boundaries and partition cell number.
170 
171 TH2Poly::TH2Poly(const char *name,const char *title,
172  Int_t nX, Double_t xlow, Double_t xup,
173  Int_t nY, Double_t ylow, Double_t yup)
174 {
175  Initialize(xlow, xup, ylow, yup, nX, nY);
176  SetName(name);
177  SetTitle(title);
178  SetFloat(kFALSE);
179 }
180 
181 ////////////////////////////////////////////////////////////////////////////////
182 /// Destructor.
183 
184 TH2Poly::~TH2Poly()
185 {
186  delete[] fCells;
187  delete[] fIsEmpty;
188  delete[] fCompletelyInside;
189  // delete at the end the bin List since it owns the objects
190  delete fBins;
191 }
192 
193 ////////////////////////////////////////////////////////////////////////////////
194 /// Create appropriate histogram bin.
195 /// e.g. TH2Poly creates TH2PolyBin,
196 /// TProfile2Poly creates TProfile2PolyBin
197 /// This is done so that TH2Poly::AddBin does not have to be duplicated,
198 /// but only create needs to be reimplemented for additional histogram types
199 
200 TH2PolyBin *TH2Poly::CreateBin(TObject *poly)
201 {
202  if (!poly) return 0;
203 
204  if (fBins == 0) {
205  fBins = new TList();
206  fBins->SetOwner();
207  }
208 
209  fNcells++;
210  Int_t ibin = fNcells - kNOverflow;
211  // if structure fsumw2 is created extend it
212  if (fSumw2.fN) fSumw2.Set(fNcells);
213  return new TH2PolyBin(poly, ibin);
214 }
215 
216 ////////////////////////////////////////////////////////////////////////////////
217 /// Adds a new bin to the histogram. It can be any object having the method
218 /// IsInside(). It returns the bin number in the histogram. It returns 0 if
219 /// it failed to add. To allow the histogram limits to expand when a bin
220 /// outside the limits is added, call SetFloat() before adding the bin.
221 
222 Int_t TH2Poly::AddBin(TObject *poly)
223 {
224  auto *bin = CreateBin(poly);
225  Int_t ibin = fNcells-kNOverflow;
226  if(!bin) return 0;
227 
228  // If the bin lies outside histogram boundaries, then extends the boundaries.
229  // Also changes the partition information accordingly
230  Bool_t flag = kFALSE;
231  if (fFloat) {
232  if (fXaxis.GetXmin() > bin->GetXMin()) {
233  fXaxis.Set(100, bin->GetXMin(), fXaxis.GetXmax());
234  flag = kTRUE;
235  }
236  if (fXaxis.GetXmax() < bin->GetXMax()) {
237  fXaxis.Set(100, fXaxis.GetXmin(), bin->GetXMax());
238  flag = kTRUE;
239  }
240  if (fYaxis.GetXmin() > bin->GetYMin()) {
241  fYaxis.Set(100, bin->GetYMin(), fYaxis.GetXmax());
242  flag = kTRUE;
243  }
244  if (fYaxis.GetXmax() < bin->GetYMax()) {
245  fYaxis.Set(100, fYaxis.GetXmin(), bin->GetYMax());
246  flag = kTRUE;
247  }
248  if (flag) ChangePartition(fCellX, fCellY);
249  } else {
250  /*Implement polygon clipping code here*/
251  }
252 
253  fBins->Add((TObject*) bin);
254  SetNewBinAdded(kTRUE);
255 
256  // Adds the bin to the partition matrix
257  AddBinToPartition(bin);
258 
259  return ibin;
260 }
261 
262 ////////////////////////////////////////////////////////////////////////////////
263 /// Adds a new bin to the histogram. The number of vertices and their (x,y)
264 /// coordinates are required as input. It returns the bin number in the
265 /// histogram.
266 
267 Int_t TH2Poly::AddBin(Int_t n, const Double_t *x, const Double_t *y)
268 {
269  TGraph *g = new TGraph(n, x, y);
270  Int_t bin = AddBin(g);
271  return bin;
272 }
273 
274 ////////////////////////////////////////////////////////////////////////////////
275 /// Add a new bin to the histogram. The bin shape is a rectangle.
276 /// It returns the bin number of the bin in the histogram.
277 
278 Int_t TH2Poly::AddBin(Double_t x1, Double_t y1, Double_t x2, Double_t y2)
279 {
280  Double_t x[] = {x1, x1, x2, x2, x1};
281  Double_t y[] = {y1, y2, y2, y1, y1};
282  TGraph *g = new TGraph(5, x, y);
283  Int_t bin = AddBin(g);
284  return bin;
285 }
286 
287 ////////////////////////////////////////////////////////////////////////////////
288 /// Performs the operation: this = this + c1*h1.
289 
290 Bool_t TH2Poly::Add(const TH1 *h1, Double_t c1)
291 {
292  Int_t bin;
293 
294  TH2Poly *h1p = (TH2Poly *)h1;
295 
296  // Check if number of bins is the same.
297  if (h1p->GetNumberOfBins() != GetNumberOfBins()) {
298  Error("Add", "Attempt to add histograms with different number of bins");
299  return kFALSE;
300  }
301 
302  // Check if the bins are the same.
303  TList *h1pBins = h1p->GetBins();
304  TH2PolyBin *thisBin, *h1pBin;
305  for (bin = 1; bin <= GetNumberOfBins(); bin++) {
306  thisBin = (TH2PolyBin *)fBins->At(bin - 1);
307  h1pBin = (TH2PolyBin *)h1pBins->At(bin - 1);
308  if (thisBin->GetXMin() != h1pBin->GetXMin() ||
309  thisBin->GetXMax() != h1pBin->GetXMax() ||
310  thisBin->GetYMin() != h1pBin->GetYMin() ||
311  thisBin->GetYMax() != h1pBin->GetYMax()) {
312  Error("Add", "Attempt to add histograms with different bin limits");
313  return kFALSE;
314  }
315  }
316 
317 
318  // Create Sumw2 if h1p has Sumw2 set
319  if (fSumw2.fN == 0 && h1p->GetSumw2N() != 0) Sumw2();
320 
321  // statistics can be preserved only in case of positive coefficients
322  // otherwise with negative c1 (histogram subtraction) one risks to get negative variances
323  Bool_t resetStats = (c1 < 0);
324  Double_t s1[kNstat] = {0};
325  Double_t s2[kNstat] = {0};
326  if (!resetStats) {
327  // need to initialize to zero s1 and s2 since
328  // GetStats fills only used elements depending on dimension and type
329  GetStats(s1);
330  h1->GetStats(s2);
331  }
332  // get number of entries now because afterwards UpdateBinContent will change it
333  Double_t entries = TMath::Abs( GetEntries() + c1 * h1->GetEntries() );
334 
335 
336  // Perform the Add.
337  Double_t factor = 1;
338  if (h1p->GetNormFactor() != 0)
339  factor = h1p->GetNormFactor() / h1p->GetSumOfWeights();
340  for (bin = 0; bin < fNcells; bin++) {
341  Double_t y = this->RetrieveBinContent(bin) + c1 * h1p->RetrieveBinContent(bin);
342  UpdateBinContent(bin, y);
343  if (fSumw2.fN) {
344  Double_t esq = factor * factor * h1p->GetBinErrorSqUnchecked(bin);
345  fSumw2.fArray[bin] += c1 * c1 * factor * factor * esq;
346  }
347  }
348 
349  // update statistics (do here to avoid changes by SetBinContent)
350  if (resetStats) {
351  // statistics need to be reset in case coefficient are negative
352  ResetStats();
353  } else {
354  for (Int_t i = 0; i < kNstat; i++) {
355  if (i == 1) s1[i] += c1 * c1 * s2[i];
356  else s1[i] += c1 * s2[i];
357  }
358  PutStats(s1);
359  SetEntries(entries);
360  }
361  return kTRUE;
362 }
363 
364 ////////////////////////////////////////////////////////////////////////////////
365 /// Adds the input bin into the partition cell matrix. This method is called
366 /// in AddBin() and ChangePartition().
367 
368 void TH2Poly::AddBinToPartition(TH2PolyBin *bin)
369 {
370  // Cell Info
371  Int_t nl, nr, mb, mt; // Max/min indices of the cells that contain the bin
372  Double_t xclipl, xclipr, yclipb, yclipt; // x and y coordinates of a cell
373  Double_t binXmax, binXmin, binYmax, binYmin; // The max/min bin coordinates
374 
375  binXmax = bin->GetXMax();
376  binXmin = bin->GetXMin();
377  binYmax = bin->GetYMax();
378  binYmin = bin->GetYMin();
379  nl = (Int_t)(floor((binXmin - fXaxis.GetXmin())/fStepX));
380  nr = (Int_t)(floor((binXmax - fXaxis.GetXmin())/fStepX));
381  mb = (Int_t)(floor((binYmin - fYaxis.GetXmin())/fStepY));
382  mt = (Int_t)(floor((binYmax - fYaxis.GetXmin())/fStepY));
383 
384  // Make sure the array indices are correct.
385  if (nr>=fCellX) nr = fCellX-1;
386  if (mt>=fCellY) mt = fCellY-1;
387  if (nl<0) nl = 0;
388  if (mb<0) mb = 0;
389 
390  // number of cells in the grid
391  //N.B. not to be confused with fNcells (the number of bins) !
392  fNCells = fCellX*fCellY;
393 
394  // Loop over all cells
395  for (int i = nl; i <= nr; i++) {
396  xclipl = fXaxis.GetXmin() + i*fStepX;
397  xclipr = xclipl + fStepX;
398  for (int j = mb; j <= mt; j++) {
399  yclipb = fYaxis.GetXmin() + j*fStepY;
400  yclipt = yclipb + fStepY;
401 
402  // If the bin is completely inside the cell,
403  // add that bin to the cell then return
404  if ((binXmin >= xclipl) && (binXmax <= xclipr) &&
405  (binYmax <= yclipt) && (binYmin >= yclipb)){
406  fCells[i + j*fCellX].Add((TObject*) bin);
407  fIsEmpty[i + j*fCellX] = kFALSE; // Makes the cell non-empty
408  return;
409  }
410 
411  // If any of the sides of the cell intersect with any side of the bin,
412  // add that bin then continue
413  if (IsIntersecting(bin, xclipl, xclipr, yclipb, yclipt)) {
414  fCells[i + j*fCellX].Add((TObject*) bin);
415  fIsEmpty[i + j*fCellX] = kFALSE; // Makes the cell non-empty
416  continue;
417  }
418  // If a corner of the cell is inside the bin and since there is no
419  // intersection, then that cell completely inside the bin.
420  if((bin->IsInside(xclipl,yclipb)) || (bin->IsInside(xclipl,yclipt))){
421  fCells[i + j*fCellX].Add((TObject*) bin);
422  fIsEmpty[i + j*fCellX] = kFALSE; // Makes the cell non-empty
423  fCompletelyInside[i + fCellX*j] = kTRUE;
424  continue;
425  }
426  if((bin->IsInside(xclipr,yclipb)) || (bin->IsInside(xclipr,yclipt))){
427  fCells[i + j*fCellX].Add((TObject*) bin);
428  fIsEmpty[i + j*fCellX] = kFALSE; // Makes the cell non-empty
429  fCompletelyInside[i + fCellX*j] = kTRUE;
430  continue;
431  }
432  }
433  }
434 }
435 
436 ////////////////////////////////////////////////////////////////////////////////
437 /// Changes the number of partition cells in the histogram.
438 /// Deletes the old partition and constructs a new one.
439 
440 void TH2Poly::ChangePartition(Int_t n, Int_t m)
441 {
442  fCellX = n; // Set the number of cells
443  fCellY = m; // Set the number of cells
444 
445  delete [] fCells; // Deletes the old partition
446 
447  // number of cells in the grid
448  //N.B. not to be confused with fNcells (the number of bins) !
449  fNCells = fCellX*fCellY;
450  fCells = new TList [fNCells]; // Sets an empty partition
451 
452  fStepX = (fXaxis.GetXmax() - fXaxis.GetXmin())/fCellX;
453  fStepY = (fYaxis.GetXmax() - fYaxis.GetXmin())/fCellY;
454 
455  delete [] fIsEmpty;
456  delete [] fCompletelyInside;
457  fIsEmpty = new Bool_t [fNCells];
458  fCompletelyInside = new Bool_t [fNCells];
459 
460  // Initializes the flags
461  for (int i = 0; i<fNCells; i++) {
462  fIsEmpty[i] = kTRUE;
463  fCompletelyInside[i] = kFALSE;
464  }
465 
466  // TList iterator
467  TIter next(fBins);
468  TObject *obj;
469 
470  while((obj = next())){ // Loop over bins and add them to the partition
471  AddBinToPartition((TH2PolyBin*) obj);
472  }
473 }
474 
475 ////////////////////////////////////////////////////////////////////////////////
476 /// Make a complete copy of the underlying object. If 'newname' is set,
477 /// the copy's name will be set to that name.
478 
479 TObject* TH2Poly::Clone(const char* newname) const
480 {
481  // TH1::Clone relies on ::Copy to implemented by the derived class.
482  // Until this is implemented, revert to the much slower default version
483  // (and possibly non-thread safe).
484 
485  return TNamed::Clone(newname);
486 }
487 
488 ////////////////////////////////////////////////////////////////////////////////
489 /// Clears the contents of all bins in the histogram.
490 
491 void TH2Poly::ClearBinContents()
492 {
493  TIter next(fBins);
494  TObject *obj;
495  TH2PolyBin *bin;
496 
497  // Clears the bin contents
498  while ((obj = next())) {
499  bin = (TH2PolyBin*) obj;
500  bin->ClearContent();
501  }
502 
503  // Clears the statistics
504  fTsumw = 0;
505  fTsumw2 = 0;
506  fTsumwx = 0;
507  fTsumwx2 = 0;
508  fTsumwy = 0;
509  fTsumwy2 = 0;
510  fEntries = 0;
511 }
512 
513 ////////////////////////////////////////////////////////////////////////////////
514 /// Reset this histogram: contents, errors, etc.
515 
516 void TH2Poly::Reset(Option_t *opt)
517 {
518  TIter next(fBins);
519  TObject *obj;
520  TH2PolyBin *bin;
521 
522  // Clears the bin contents
523  while ((obj = next())) {
524  bin = (TH2PolyBin*) obj;
525  bin->ClearContent();
526  }
527 
528  TH2::Reset(opt);
529 }
530 
531 ////////////////////////////////////////////////////////////////////////////////
532 /// Returns the bin number of the bin at the given coordinate. -1 to -9 are
533 /// the overflow and underflow bins. overflow bin -5 is the unbinned areas in
534 /// the histogram (also called "the sea"). The third parameter can be left
535 /// blank.
536 /// The overflow/underflow bins are:
537 ///~~~ {.cpp}
538 /// -1 | -2 | -3
539 /// -------------
540 /// -4 | -5 | -6
541 /// -------------
542 /// -7 | -8 | -9
543 ///~~~
544 /// where -5 means is the "sea" bin (i.e. unbinned areas)
545 
546 Int_t TH2Poly::FindBin(Double_t x, Double_t y, Double_t)
547 {
548 
549  // Checks for overflow/underflow
550  Int_t overflow = 0;
551  if (y > fYaxis.GetXmax()) overflow += -1;
552  else if (y > fYaxis.GetXmin()) overflow += -4;
553  else overflow += -7;
554  if (x > fXaxis.GetXmax()) overflow += -2;
555  else if (x > fXaxis.GetXmin()) overflow += -1;
556  if (overflow != -5) return overflow;
557 
558  // Finds the cell (x,y) coordinates belong to
559  Int_t n = (Int_t)(floor((x-fXaxis.GetXmin())/fStepX));
560  Int_t m = (Int_t)(floor((y-fYaxis.GetXmin())/fStepY));
561 
562  // Make sure the array indices are correct.
563  if (n>=fCellX) n = fCellX-1;
564  if (m>=fCellY) m = fCellY-1;
565  if (n<0) n = 0;
566  if (m<0) m = 0;
567 
568  if (fIsEmpty[n+fCellX*m]) return -5;
569 
570  TH2PolyBin *bin;
571 
572  TIter next(&fCells[n+fCellX*m]);
573  TObject *obj;
574 
575  // Search for the bin in the cell
576  while ((obj=next())) {
577  bin = (TH2PolyBin*)obj;
578  if (bin->IsInside(x,y)) return bin->GetBinNumber();
579  }
580 
581  // If the search has not returned a bin, the point must be on "the sea"
582  return -5;
583 }
584 
585 ////////////////////////////////////////////////////////////////////////////////
586 /// Increment the bin containing (x,y) by 1.
587 /// Uses the partitioning algorithm.
588 
589 Int_t TH2Poly::Fill(Double_t x, Double_t y)
590 {
591  return Fill(x, y, 1.0);
592 }
593 
594 ////////////////////////////////////////////////////////////////////////////////
595 /// Increment the bin containing (x,y) by w.
596 /// Uses the partitioning algorithm.
597 
598 Int_t TH2Poly::Fill(Double_t x, Double_t y, Double_t w)
599 {
600  // see GetBinCOntent for definition of overflow bins
601  // in case of weighted events store weight square in fSumw2.fArray
602  // but with this indexing:
603  // fSumw2.fArray[0:kNOverflow-1] : sum of weight squares for the overflow bins
604  // fSumw2.fArray[kNOverflow:fNcells] : sum of weight squares for the standard bins
605  // where fNcells = kNOverflow + Number of bins. kNOverflow=9
606 
607  if (fNcells <= kNOverflow) return 0;
608 
609  // create sum of weight square array if weights are different than 1
610  if (!fSumw2.fN && w != 1.0 && !TestBit(TH1::kIsNotW) ) Sumw2();
611 
612  Int_t overflow = 0;
613  if (y > fYaxis.GetXmax()) overflow += -1;
614  else if (y > fYaxis.GetXmin()) overflow += -4;
615  else overflow += -7;
616  if (x > fXaxis.GetXmax()) overflow += -2;
617  else if(x > fXaxis.GetXmin()) overflow += -1;
618  if (overflow != -5) {
619  fOverflow[-overflow - 1]+= w;
620  if (fSumw2.fN) fSumw2.fArray[-overflow - 1] += w*w;
621  return overflow;
622  }
623 
624  // Finds the cell (x,y) coordinates belong to
625  Int_t n = (Int_t)(floor((x-fXaxis.GetXmin())/fStepX));
626  Int_t m = (Int_t)(floor((y-fYaxis.GetXmin())/fStepY));
627 
628  // Make sure the array indices are correct.
629  if (n>=fCellX) n = fCellX-1;
630  if (m>=fCellY) m = fCellY-1;
631  if (n<0) n = 0;
632  if (m<0) m = 0;
633 
634  if (fIsEmpty[n+fCellX*m]) {
635  fOverflow[4]+= w;
636  if (fSumw2.fN) fSumw2.fArray[4] += w*w;
637  return -5;
638  }
639 
640  TH2PolyBin *bin;
641  Int_t bi;
642 
643  TIter next(&fCells[n+fCellX*m]);
644  TObject *obj;
645 
646  while ((obj=next())) {
647  bin = (TH2PolyBin*)obj;
648  // needs to account offset in array for overflow bins
649  bi = bin->GetBinNumber()-1+kNOverflow;
650  if (bin->IsInside(x,y)) {
651  bin->Fill(w);
652 
653  // Statistics
654  fTsumw = fTsumw + w;
655  fTsumw2 = fTsumw2 + w*w;
656  fTsumwx = fTsumwx + w*x;
657  fTsumwx2 = fTsumwx2 + w*x*x;
658  fTsumwy = fTsumwy + w*y;
659  fTsumwy2 = fTsumwy2 + w*y*y;
660  if (fSumw2.fN) {
661  assert(bi < fSumw2.fN);
662  fSumw2.fArray[bi] += w*w;
663  }
664  fEntries++;
665 
666  SetBinContentChanged(kTRUE);
667 
668  return bin->GetBinNumber();
669  }
670  }
671 
672  fOverflow[4]+= w;
673  if (fSumw2.fN) fSumw2.fArray[4] += w*w;
674  return -5;
675 }
676 
677 ////////////////////////////////////////////////////////////////////////////////
678 /// Increment the bin named "name" by w.
679 
680 Int_t TH2Poly::Fill(const char* name, Double_t w)
681 {
682  TString sname(name);
683 
684  TIter next(fBins);
685  TObject *obj;
686  TH2PolyBin *bin;
687 
688  while ((obj = next())) {
689  bin = (TH2PolyBin*) obj;
690  if (!sname.CompareTo(bin->GetPolygon()->GetName())) {
691  bin->Fill(w);
692  fEntries++;
693  SetBinContentChanged(kTRUE);
694  return bin->GetBinNumber();
695  }
696  }
697 
698  return 0;
699 }
700 
701 ////////////////////////////////////////////////////////////////////////////////
702 /// Fills a 2-D histogram with an array of values and weights.
703 ///
704 /// \param [in] ntimes: number of entries in arrays x and w
705 /// (array size must be ntimes*stride)
706 /// \param [in] x: array of x values to be histogrammed
707 /// \param [in] y: array of y values to be histogrammed
708 /// \param [in] w: array of weights
709 /// \param [in] stride: step size through arrays x, y and w
710 
711 void TH2Poly::FillN(Int_t ntimes, const Double_t* x, const Double_t* y,
712  const Double_t* w, Int_t stride)
713 {
714  for (int i = 0; i < ntimes; i += stride) {
715  Fill(x[i], y[i], w[i]);
716  }
717 }
718 
719 ////////////////////////////////////////////////////////////////////////////////
720 /// Returns the integral of bin contents.
721 /// By default the integral is computed as the sum of bin contents.
722 /// If option "width" or "area" is specified, the integral is the sum of
723 /// the bin contents multiplied by the area of the bin.
724 
725 Double_t TH2Poly::Integral(Option_t* option) const
726 {
727  TString opt = option;
728  opt.ToLower();
729 
730  if ((opt.Contains("width")) || (opt.Contains("area"))) {
731  Double_t w;
732  Double_t integral = 0.;
733 
734  TIter next(fBins);
735  TObject *obj;
736  TH2PolyBin *bin;
737  while ((obj=next())) {
738  bin = (TH2PolyBin*) obj;
739  w = bin->GetArea();
740  integral += w*(bin->GetContent());
741  }
742 
743  return integral;
744  } else {
745  return fTsumw;
746  }
747 }
748 
749 ////////////////////////////////////////////////////////////////////////////////
750 /// Returns the content of the input bin
751 /// For the overflow/underflow/sea bins:
752 ///~~~ {.cpp}
753 /// -1 | -2 | -3
754 /// ---+----+----
755 /// -4 | -5 | -6
756 /// ---+----+----
757 /// -7 | -8 | -9
758 ///~~~
759 /// where -5 is the "sea" bin (i.e. unbinned areas)
760 
761 Double_t TH2Poly::GetBinContent(Int_t bin) const
762 {
763  if (bin > GetNumberOfBins() || bin == 0 || bin < -kNOverflow) return 0;
764  if (bin<0) return fOverflow[-bin - 1];
765  return ((TH2PolyBin*) fBins->At(bin-1))->GetContent();
766 }
767 
768 ////////////////////////////////////////////////////////////////////////////////
769 /// Returns the value of error associated to bin number bin.
770 /// If the sum of squares of weights has been defined (via Sumw2),
771 /// this function returns the sqrt(sum of w2).
772 /// otherwise it returns the sqrt(contents) for this bin.
773 /// Bins are in range [1:nbins] and for bin < 0 in range [-9:-1] it returns errors for overflow bins.
774 /// See also TH2Poly::GetBinContent
775 
776 Double_t TH2Poly::GetBinError(Int_t bin) const
777 {
778  if (bin == 0 || bin > GetNumberOfBins() || bin < - kNOverflow) return 0;
779  if (fBuffer) ((TH1*)this)->BufferEmpty();
780  // in case of weighted events the sum of the weights are stored in a different way than
781  // a normal histogram
782  // fSumw2.fArray[0:kNOverflow-1] : sum of weight squares for the overflow bins (
783  // fSumw2.fArray[kNOverflow:fNcells] : sum of weight squares for the standard bins
784  // fNcells = kNOverflow (9) + Number of bins
785  if (fSumw2.fN) {
786  Int_t binIndex = (bin > 0) ? bin+kNOverflow-1 : -(bin+1);
787  Double_t err2 = fSumw2.fArray[binIndex];
788  return TMath::Sqrt(err2);
789  }
790  Double_t error2 = TMath::Abs(GetBinContent(bin));
791  return TMath::Sqrt(error2);
792 }
793 
794 ////////////////////////////////////////////////////////////////////////////////
795 /// Set the bin Error.
796 /// Re-implementation for TH2Poly given the different bin indexing in the
797 /// stored squared error array.
798 /// See also notes in TH1::SetBinError
799 ///
800 /// Bins are in range [1:nbins] and for bin < 0 in the range [-9:-1] the errors is set for the overflow bins
801 
802 
803 void TH2Poly::SetBinError(Int_t bin, Double_t error)
804 {
805  if (bin == 0 || bin > GetNumberOfBins() || bin < - kNOverflow) return;
806  if (!fSumw2.fN) Sumw2();
807  SetBinErrorOption(kNormal);
808  // see comment in GetBinError for special convention of bin index in fSumw2 array
809  Int_t binIndex = (bin > 0) ? bin+kNOverflow-1 : -(bin+1);
810  fSumw2.fArray[binIndex] = error * error;
811 }
812 
813 
814 
815 ////////////////////////////////////////////////////////////////////////////////
816 /// Returns the bin name.
817 
818 const char *TH2Poly::GetBinName(Int_t bin) const
819 {
820  if (bin > GetNumberOfBins()) return "";
821  if (bin < 0) return "";
822  return ((TH2PolyBin*) fBins->At(bin-1))->GetPolygon()->GetName();
823 }
824 
825 ////////////////////////////////////////////////////////////////////////////////
826 /// Returns the bin title.
827 
828 const char *TH2Poly::GetBinTitle(Int_t bin) const
829 {
830  if (bin > GetNumberOfBins()) return "";
831  if (bin < 0) return "";
832  return ((TH2PolyBin*) fBins->At(bin-1))->GetPolygon()->GetTitle();
833 }
834 
835 ////////////////////////////////////////////////////////////////////////////////
836 /// Returns the maximum value of the histogram.
837 
838 Double_t TH2Poly::GetMaximum() const
839 {
840  if (fNcells <= kNOverflow) return 0;
841  if (fMaximum != -1111) return fMaximum;
842 
843  TH2PolyBin *b;
844 
845  TIter next(fBins);
846  TObject *obj;
847  Double_t max,c;
848 
849  max = ((TH2PolyBin*) next())->GetContent();
850 
851  while ((obj=next())) {
852  b = (TH2PolyBin*)obj;
853  c = b->GetContent();
854  if (c>max) max = c;
855  }
856  return max;
857 }
858 
859 ////////////////////////////////////////////////////////////////////////////////
860 /// Returns the maximum value of the histogram that is less than maxval.
861 
862 Double_t TH2Poly::GetMaximum(Double_t maxval) const
863 {
864  if (fNcells <= kNOverflow) return 0;
865  if (fMaximum != -1111) return fMaximum;
866 
867  TH2PolyBin *b;
868 
869  TIter next(fBins);
870  TObject *obj;
871  Double_t max,c;
872 
873  max = ((TH2PolyBin*) next())->GetContent();
874 
875  while ((obj=next())) {
876  b = (TH2PolyBin*)obj;
877  c = b->GetContent();
878  if (c>max && c<maxval) max=c;
879  }
880  return max;
881 }
882 
883 ////////////////////////////////////////////////////////////////////////////////
884 /// Returns the minimum value of the histogram.
885 
886 Double_t TH2Poly::GetMinimum() const
887 {
888  if (fNcells <= kNOverflow) return 0;
889  if (fMinimum != -1111) return fMinimum;
890 
891  TH2PolyBin *b;
892 
893  TIter next(fBins);
894  TObject *obj;
895  Double_t min,c;
896 
897  min = ((TH2PolyBin*) next())->GetContent();
898 
899  while ((obj=next())) {
900  b = (TH2PolyBin*)obj;
901  c = b->GetContent();
902  if (c<min) min=c;
903  }
904  return min;
905 }
906 
907 ////////////////////////////////////////////////////////////////////////////////
908 /// Returns the minimum value of the histogram that is greater than minval.
909 
910 Double_t TH2Poly::GetMinimum(Double_t minval) const
911 {
912  if (fNcells <= kNOverflow) return 0;
913  if (fMinimum != -1111) return fMinimum;
914 
915  TH2PolyBin *b;
916 
917  TIter next(fBins);
918  TObject *obj;
919  Double_t min,c;
920 
921  min = ((TH2PolyBin*) next())->GetContent();
922 
923  while ((obj=next())) {
924  b = (TH2PolyBin*)obj;
925  c = b->GetContent();
926  if (c<min && c>minval) min=c;
927  }
928  return min;
929 }
930 
931 ////////////////////////////////////////////////////////////////////////////////
932 /// Bins the histogram using a honeycomb structure
933 
934 void TH2Poly::Honeycomb(Double_t xstart, Double_t ystart, Double_t a,
935  Int_t k, Int_t s)
936 {
937  // Add the bins
938  Double_t numberOfHexagonsInTheRow;
939  Double_t x[6], y[6];
940  Double_t xloop, yloop, xtemp;
941  xloop = xstart; yloop = ystart + a/2.0;
942  for (int sCounter = 0; sCounter < s; sCounter++) {
943 
944  xtemp = xloop; // Resets the temp variable
945 
946  // Determine the number of hexagons in that row
947  if(sCounter%2 == 0){numberOfHexagonsInTheRow = k;}
948  else{numberOfHexagonsInTheRow = k - 1;}
949 
950  for (int kCounter = 0; kCounter < numberOfHexagonsInTheRow; kCounter++) {
951 
952  // Go around the hexagon
953  x[0] = xtemp;
954  y[0] = yloop;
955  x[1] = x[0];
956  y[1] = y[0] + a;
957  x[2] = x[1] + a*TMath::Sqrt(3)/2.0;
958  y[2] = y[1] + a/2.0;
959  x[3] = x[2] + a*TMath::Sqrt(3)/2.0;
960  y[3] = y[1];
961  x[4] = x[3];
962  y[4] = y[0];
963  x[5] = x[2];
964  y[5] = y[4] - a/2.0;
965 
966  this->AddBin(6, x, y);
967 
968  // Go right
969  xtemp += a*TMath::Sqrt(3);
970  }
971 
972  // Increment the starting position
973  if (sCounter%2 == 0) xloop += a*TMath::Sqrt(3)/2.0;
974  else xloop -= a*TMath::Sqrt(3)/2.0;
975  yloop += 1.5*a;
976  }
977 }
978 
979 ////////////////////////////////////////////////////////////////////////////////
980 /// Initializes the TH2Poly object. This method is called by the constructor.
981 
982 void TH2Poly::Initialize(Double_t xlow, Double_t xup,
983  Double_t ylow, Double_t yup, Int_t n, Int_t m)
984 {
985  Int_t i;
986  fDimension = 2; //The dimension of the histogram
987 
988  fBins = 0;
989  fNcells = kNOverflow;
990 
991  // Sets the boundaries of the histogram
992  fXaxis.Set(100, xlow, xup);
993  fYaxis.Set(100, ylow, yup);
994 
995  for (i=0; i<9; i++) fOverflow[i] = 0.;
996 
997  // Statistics
998  fEntries = 0; // The total number of entries
999  fTsumw = 0.; // Total amount of content in the histogram
1000  fTsumw2 = 0.; // Sum square of the weights
1001  fTsumwx = 0.; // Weighted sum of x coordinates
1002  fTsumwx2 = 0.; // Weighted sum of the squares of x coordinates
1003  fTsumwy2 = 0.; // Weighted sum of the squares of y coordinates
1004  fTsumwy = 0.; // Weighted sum of y coordinates
1005 
1006  fCellX = n; // Set the number of cells to default
1007  fCellY = m; // Set the number of cells to default
1008 
1009  // number of cells in the grid
1010  //N.B. not to be confused with fNcells (the number of bins) !
1011  fNCells = fCellX*fCellY;
1012  fCells = new TList [fNCells]; // Sets an empty partition
1013  fStepX = (fXaxis.GetXmax() - fXaxis.GetXmin())/fCellX; // Cell width
1014  fStepY = (fYaxis.GetXmax() - fYaxis.GetXmin())/fCellY; // Cell height
1015 
1016  fIsEmpty = new Bool_t [fNCells]; // Empty partition
1017  fCompletelyInside = new Bool_t [fNCells]; // Cell is completely inside bin
1018 
1019  for (i = 0; i<fNCells; i++) { // Initializes the flags
1020  fIsEmpty[i] = kTRUE;
1021  fCompletelyInside[i] = kFALSE;
1022  }
1023 
1024  // 3D Painter flags
1025  SetNewBinAdded(kFALSE);
1026  SetBinContentChanged(kFALSE);
1027 }
1028 
1029 ////////////////////////////////////////////////////////////////////////////////
1030 /// Returns kTRUE if the input bin is intersecting with the
1031 /// input rectangle (xclipl, xclipr, yclipb, yclipt)
1032 
1033 Bool_t TH2Poly::IsIntersecting(TH2PolyBin *bin,
1034  Double_t xclipl, Double_t xclipr,
1035  Double_t yclipb, Double_t yclipt)
1036 {
1037  Int_t gn;
1038  Double_t *gx;
1039  Double_t *gy;
1040  Bool_t inter = kFALSE;
1041  TObject *poly = bin->GetPolygon();
1042 
1043  if (poly->IsA() == TGraph::Class()) {
1044  TGraph *g = (TGraph*)poly;
1045  gx = g->GetX();
1046  gy = g->GetY();
1047  gn = g->GetN();
1048  inter = IsIntersectingPolygon(gn, gx, gy, xclipl, xclipr, yclipb, yclipt);
1049  }
1050 
1051  if (poly->IsA() == TMultiGraph::Class()) {
1052  TMultiGraph *mg = (TMultiGraph*)poly;
1053  TList *gl = mg->GetListOfGraphs();
1054  if (!gl) return inter;
1055  TGraph *g;
1056  TIter next(gl);
1057  while ((g = (TGraph*) next())) {
1058  gx = g->GetX();
1059  gy = g->GetY();
1060  gn = g->GetN();
1061  inter = IsIntersectingPolygon(gn, gx, gy, xclipl, xclipr,
1062  yclipb, yclipt);
1063  if (inter) return inter;
1064  }
1065  }
1066 
1067  return inter;
1068 }
1069 
1070 ////////////////////////////////////////////////////////////////////////////////
1071 /// Returns kTRUE if the input polygon (bn, x, y) is intersecting with the
1072 /// input rectangle (xclipl, xclipr, yclipb, yclipt)
1073 
1074 Bool_t TH2Poly::IsIntersectingPolygon(Int_t bn, Double_t *x, Double_t *y,
1075  Double_t xclipl, Double_t xclipr,
1076  Double_t yclipb, Double_t yclipt)
1077 {
1078  Bool_t p0R, p0L, p0T, p0B, p0xM, p0yM, p1R, p1L, p1T;
1079  Bool_t p1B, p1xM, p1yM, p0In, p1In;
1080 
1081  for (int counter = 0; counter < (bn-1); counter++) {
1082  // If both are on the same side, return kFALSE
1083  p0L = x[counter] <= xclipl; // Point 0 is on the left
1084  p1L = x[counter + 1] <= xclipl; // Point 1 is on the left
1085  if (p0L && p1L) continue;
1086  p0R = x[counter] >= xclipr; // Point 0 is on the right
1087  p1R = x[counter + 1] >= xclipr; // Point 1 is on the right
1088  if (p0R && p1R) continue;
1089  p0T = y[counter] >= yclipt; // Point 0 is at the top
1090  p1T = y[counter + 1] >= yclipt; // Point 1 is at the top
1091  if (p0T && p1T) continue;
1092  p0B = y[counter] <= yclipb; // Point 0 is at the bottom
1093  p1B = y[counter + 1] <= yclipb; // Point 1 is at the bottom
1094  if (p0B && p1B) continue;
1095 
1096  // Checks to see if any are inside
1097  p0xM = !p0R && !p0L; // Point 0 is inside along x
1098  p0yM = !p0T && !p0B; // Point 1 is inside along x
1099  p1xM = !p1R && !p1L; // Point 0 is inside along y
1100  p1yM = !p1T && !p1B; // Point 1 is inside along y
1101  p0In = p0xM && p0yM; // Point 0 is inside
1102  p1In = p1xM && p1yM; // Point 1 is inside
1103  if (p0In) {
1104  if (p1In) continue;
1105  return kTRUE;
1106  } else {
1107  if (p1In) return kTRUE;
1108  }
1109 
1110  // We know by now that the points are not in the same side and not inside.
1111 
1112  // Checks to see if they are opposite
1113 
1114  if (p0xM && p1xM) return kTRUE;
1115  if (p0yM && p1yM) return kTRUE;
1116 
1117  // We now know that the points are in different x and y indices
1118 
1119  Double_t xcoord[3], ycoord[3];
1120  xcoord[0] = x[counter];
1121  xcoord[1] = x[counter + 1];
1122  ycoord[0] = y[counter];
1123  ycoord[1] = y[counter + 1];
1124 
1125  if (p0L) {
1126  if(p1T){
1127  xcoord[2] = xclipl;
1128  ycoord[2] = yclipb;
1129  if((TMath::IsInside(xclipl, yclipt, 3, xcoord, ycoord)) ||
1130  (TMath::IsInside(xclipr, yclipb, 3, xcoord, ycoord))) continue;
1131  else return kTRUE;
1132  } else if (p1B) {
1133  xcoord[2] = xclipl;
1134  ycoord[2] = yclipt;
1135  if((TMath::IsInside(xclipl, yclipb, 3, xcoord, ycoord)) ||
1136  (TMath::IsInside(xclipr, yclipt, 3, xcoord, ycoord))) continue;
1137  else return kTRUE;
1138  } else { // p1yM
1139  if (p0T) {
1140  xcoord[2] = xclipl;
1141  ycoord[2] = yclipb;
1142  if (TMath::IsInside(xclipr, yclipt, 3, xcoord, ycoord)) continue;
1143  else return kTRUE;
1144  }
1145  if (p0B) {
1146  xcoord[2] = xclipl;
1147  ycoord[2] = yclipt;
1148  if (TMath::IsInside(xclipr, yclipb, 3, xcoord, ycoord)) continue;
1149  else return kTRUE;
1150  }
1151  }
1152  } else if (p0R) {
1153  if (p1T) {
1154  xcoord[2] = xclipl;
1155  ycoord[2] = yclipb;
1156  if ((TMath::IsInside(xclipr, yclipb, 3, xcoord, ycoord)) ||
1157  (TMath::IsInside(xclipl, yclipt, 3, xcoord, ycoord))) continue;
1158  else return kTRUE;
1159  } else if (p1B) {
1160  xcoord[2] = xclipl;
1161  ycoord[2] = yclipt;
1162  if ((TMath::IsInside(xclipl, yclipb, 3, xcoord, ycoord)) ||
1163  (TMath::IsInside(xclipr, yclipt, 3, xcoord, ycoord))) continue;
1164  else return kTRUE;
1165  } else{ // p1yM
1166  if (p0T) {
1167  xcoord[2] = xclipr;
1168  ycoord[2] = yclipb;
1169  if (TMath::IsInside(xclipl, yclipt, 3, xcoord, ycoord)) continue;
1170  else return kTRUE;
1171  }
1172  if (p0B) {
1173  xcoord[2] = xclipr;
1174  ycoord[2] = yclipt;
1175  if (TMath::IsInside(xclipl, yclipb, 3, xcoord, ycoord)) continue;
1176  else return kTRUE;
1177  }
1178  }
1179  }
1180  }
1181  return kFALSE;
1182 }
1183 
1184 ////////////////////////////////////////////////////////////////////////////////
1185 /// Merge TH2Polys
1186 /// Given the special nature of the TH2Poly, the merge is implemented in
1187 /// terms of subsequent TH2Poly::Add calls.
1188 Long64_t TH2Poly::Merge(TCollection *coll)
1189 {
1190  for (auto h2pAsObj : *coll) {
1191  if (!Add((TH1*)h2pAsObj, 1.)) {
1192  Warning("Merge", "An issue was encountered during the merge operation.");
1193  return 0L;
1194  }
1195  }
1196  return GetEntries();
1197 }
1198 
1199 ////////////////////////////////////////////////////////////////////////////////
1200 /// Save primitive as a C++ statement(s) on output stream out
1201 
1202 void TH2Poly::SavePrimitive(std::ostream &out, Option_t *option)
1203 {
1204  out <<" "<<std::endl;
1205  out <<" "<< ClassName() <<" *";
1206 
1207  //histogram pointer has by default the histogram name.
1208  //however, in case histogram has no directory, it is safer to add a
1209  //incremental suffix
1210  static Int_t hcounter = 0;
1211  TString histName = GetName();
1212  if (!fDirectory && !histName.Contains("Graph")) {
1213  hcounter++;
1214  histName += "__";
1215  histName += hcounter;
1216  }
1217  const char *hname = histName.Data();
1218 
1219  //Construct the class initialization
1220  out << hname << " = new " << ClassName() << "(\"" << hname << "\", \""
1221  << GetTitle() << "\", " << fCellX << ", " << fXaxis.GetXmin()
1222  << ", " << fXaxis.GetXmax()
1223  << ", " << fCellY << ", " << fYaxis.GetXmin() << ", "
1224  << fYaxis.GetXmax() << ");" << std::endl;
1225 
1226  // Save Bins
1227  TIter next(fBins);
1228  TObject *obj;
1229  TH2PolyBin *th2pBin;
1230 
1231  while((obj = next())){
1232  th2pBin = (TH2PolyBin*) obj;
1233  th2pBin->GetPolygon()->SavePrimitive(out,
1234  TString::Format("th2poly%s",histName.Data()));
1235  }
1236 
1237  // save bin contents
1238  out<<" "<<std::endl;
1239  Int_t bin;
1240  for (bin=1;bin<=GetNumberOfBins();bin++) {
1241  Double_t bc = GetBinContent(bin);
1242  if (bc) {
1243  out<<" "<<hname<<"->SetBinContent("<<bin<<","<<bc<<");"<<std::endl;
1244  }
1245  }
1246 
1247  // save bin errors
1248  if (fSumw2.fN) {
1249  for (bin=1;bin<=GetNumberOfBins();bin++) {
1250  Double_t be = GetBinError(bin);
1251  if (be) {
1252  out<<" "<<hname<<"->SetBinError("<<bin<<","<<be<<");"<<std::endl;
1253  }
1254  }
1255  }
1256  TH1::SavePrimitiveHelp(out, hname, option);
1257 }
1258 
1259 ////////////////////////////////////////////////////////////////////////////////
1260 /// Multiply this histogram by a constant c1.
1261 
1262 void TH2Poly::Scale(Double_t c1, Option_t*)
1263 {
1264  for( int i = 0; i < this->GetNumberOfBins(); i++ ) {
1265  this->SetBinContent(i+1, c1*this->GetBinContent(i+1));
1266  }
1267  for( int i = 0; i < kNOverflow; i++ ) {
1268  this->SetBinContent(-i-1, c1*this->GetBinContent(-i-1) );
1269  }
1270 }
1271 
1272 ////////////////////////////////////////////////////////////////////////////////
1273 /// Sets the contents of the input bin to the input content
1274 /// Negative values between -1 and -9 are for the overflows and the sea
1275 
1276 void TH2Poly::SetBinContent(Int_t bin, Double_t content)
1277 {
1278  if (bin > GetNumberOfBins() || bin == 0 || bin < -9 ) return;
1279  if (bin > 0) {
1280  ((TH2PolyBin*) fBins->At(bin-1))->SetContent(content);
1281  }
1282  else
1283  fOverflow[-bin - 1] = content;
1284 
1285  SetBinContentChanged(kTRUE);
1286  fEntries++;
1287 }
1288 
1289 ////////////////////////////////////////////////////////////////////////////////
1290 /// When set to kTRUE, allows the histogram to expand if a bin outside the
1291 /// limits is added.
1292 
1293 void TH2Poly::SetFloat(Bool_t flag)
1294 {
1295  fFloat = flag;
1296 }
1297 
1298 ////////////////////////////////////////////////////////////////////////////////
1299 /// Return "true" if the point (x,y) is inside the bin of binnr.
1300 
1301 Bool_t TH2Poly::IsInsideBin(Int_t binnr, Double_t x, Double_t y)
1302 {
1303  if (!fBins) return false;
1304  TH2PolyBin* bin = (TH2PolyBin*)fBins->At(binnr);
1305  if (!bin) return false;
1306  return bin->IsInside(x,y);
1307 }
1308 
1309 void TH2Poly::GetStats(Double_t *stats) const
1310 {
1311  stats[0] = fTsumw;
1312  stats[1] = fTsumw2;
1313  stats[2] = fTsumwx;
1314  stats[3] = fTsumwx2;
1315  stats[4] = fTsumwy;
1316  stats[5] = fTsumwy2;
1317  stats[6] = fTsumwxy;
1318 }
1319 
1320 /** \class TH2PolyBin
1321  \ingroup Hist
1322 Helper class to represent a bin in the TH2Poly histogram
1323 */
1324 
1325 ////////////////////////////////////////////////////////////////////////////////
1326 /// Default constructor.
1327 
1328 TH2PolyBin::TH2PolyBin()
1329 {
1330  fPoly = 0;
1331  fContent = 0.;
1332  fNumber = 0;
1333  fXmax = -1111;
1334  fXmin = -1111;
1335  fYmax = -1111;
1336  fYmin = -1111;
1337  fArea = 0;
1338  SetChanged(kTRUE);
1339 }
1340 
1341 ////////////////////////////////////////////////////////////////////////////////
1342 /// Normal constructor.
1343 
1344 TH2PolyBin::TH2PolyBin(TObject *poly, Int_t bin_number)
1345 {
1346  fContent = 0.;
1347  fNumber = bin_number;
1348  fArea = 0.;
1349  fPoly = poly;
1350  fXmax = -1111;
1351  fXmin = -1111;
1352  fYmax = -1111;
1353  fYmin = -1111;
1354  SetChanged(kTRUE);
1355 }
1356 
1357 ////////////////////////////////////////////////////////////////////////////////
1358 /// Destructor.
1359 
1360 TH2PolyBin::~TH2PolyBin()
1361 {
1362  if (fPoly) delete fPoly;
1363 }
1364 
1365 ////////////////////////////////////////////////////////////////////////////////
1366 /// Returns the area of the bin.
1367 
1368 Double_t TH2PolyBin::GetArea()
1369 {
1370  Int_t bn;
1371 
1372  if (fArea == 0) {
1373  if (fPoly->IsA() == TGraph::Class()) {
1374  TGraph *g = (TGraph*)fPoly;
1375  bn = g->GetN();
1376  fArea = g->Integral(0,bn-1);
1377  }
1378 
1379  if (fPoly->IsA() == TMultiGraph::Class()) {
1380  TMultiGraph *mg = (TMultiGraph*)fPoly;
1381  TList *gl = mg->GetListOfGraphs();
1382  if (!gl) return fArea;
1383  TGraph *g;
1384  TIter next(gl);
1385  while ((g = (TGraph*) next())) {
1386  bn = g->GetN();
1387  fArea = fArea + g->Integral(0,bn-1);
1388  }
1389  }
1390  }
1391 
1392  return fArea;
1393 }
1394 
1395 ////////////////////////////////////////////////////////////////////////////////
1396 /// Returns the maximum value for the x coordinates of the bin.
1397 
1398 Double_t TH2PolyBin::GetXMax()
1399 {
1400  if (fXmax != -1111) return fXmax;
1401 
1402  Int_t bn,i;
1403  Double_t *bx;
1404 
1405  if (fPoly->IsA() == TGraph::Class()) {
1406  TGraph *g = (TGraph*)fPoly;
1407  bx = g->GetX();
1408  bn = g->GetN();
1409  fXmax = bx[0];
1410  for (i=1; i<bn; i++) {if (fXmax < bx[i]) fXmax = bx[i];}
1411  }
1412 
1413  if (fPoly->IsA() == TMultiGraph::Class()) {
1414  TMultiGraph *mg = (TMultiGraph*)fPoly;
1415  TList *gl = mg->GetListOfGraphs();
1416  if (!gl) return fXmax;
1417  TGraph *g;
1418  TIter next(gl);
1419  Bool_t first = kTRUE;
1420  while ((g = (TGraph*) next())) {
1421  bx = g->GetX();
1422  bn = g->GetN();
1423  if (first) {fXmax = bx[0]; first = kFALSE;}
1424  for (i=0; i<bn; i++) {if (fXmax < bx[i]) fXmax = bx[i];}
1425  }
1426  }
1427 
1428  return fXmax;
1429 }
1430 
1431 ////////////////////////////////////////////////////////////////////////////////
1432 /// Returns the minimum value for the x coordinates of the bin.
1433 
1434 Double_t TH2PolyBin::GetXMin()
1435 {
1436  if (fXmin != -1111) return fXmin;
1437 
1438  Int_t bn,i;
1439  Double_t *bx;
1440 
1441  if (fPoly->IsA() == TGraph::Class()) {
1442  TGraph *g = (TGraph*)fPoly;
1443  bx = g->GetX();
1444  bn = g->GetN();
1445  fXmin = bx[0];
1446  for (i=1; i<bn; i++) {if (fXmin > bx[i]) fXmin = bx[i];}
1447  }
1448 
1449  if (fPoly->IsA() == TMultiGraph::Class()) {
1450  TMultiGraph *mg = (TMultiGraph*)fPoly;
1451  TList *gl = mg->GetListOfGraphs();
1452  if (!gl) return fXmin;
1453  TGraph *g;
1454  TIter next(gl);
1455  Bool_t first = kTRUE;
1456  while ((g = (TGraph*) next())) {
1457  bx = g->GetX();
1458  bn = g->GetN();
1459  if (first) {fXmin = bx[0]; first = kFALSE;}
1460  for (i=0; i<bn; i++) {if (fXmin > bx[i]) fXmin = bx[i];}
1461  }
1462  }
1463 
1464  return fXmin;
1465 }
1466 
1467 ////////////////////////////////////////////////////////////////////////////////
1468 /// Returns the maximum value for the y coordinates of the bin.
1469 
1470 Double_t TH2PolyBin::GetYMax()
1471 {
1472  if (fYmax != -1111) return fYmax;
1473 
1474  Int_t bn,i;
1475  Double_t *by;
1476 
1477  if (fPoly->IsA() == TGraph::Class()) {
1478  TGraph *g = (TGraph*)fPoly;
1479  by = g->GetY();
1480  bn = g->GetN();
1481  fYmax = by[0];
1482  for (i=1; i<bn; i++) {if (fYmax < by[i]) fYmax = by[i];}
1483  }
1484 
1485  if (fPoly->IsA() == TMultiGraph::Class()) {
1486  TMultiGraph *mg = (TMultiGraph*)fPoly;
1487  TList *gl = mg->GetListOfGraphs();
1488  if (!gl) return fYmax;
1489  TGraph *g;
1490  TIter next(gl);
1491  Bool_t first = kTRUE;
1492  while ((g = (TGraph*) next())) {
1493  by = g->GetY();
1494  bn = g->GetN();
1495  if (first) {fYmax = by[0]; first = kFALSE;}
1496  for (i=0; i<bn; i++) {if (fYmax < by[i]) fYmax = by[i];}
1497  }
1498  }
1499 
1500  return fYmax;
1501 }
1502 
1503 ////////////////////////////////////////////////////////////////////////////////
1504 /// Returns the minimum value for the y coordinates of the bin.
1505 
1506 Double_t TH2PolyBin::GetYMin()
1507 {
1508  if (fYmin != -1111) return fYmin;
1509 
1510  Int_t bn,i;
1511  Double_t *by;
1512 
1513  if (fPoly->IsA() == TGraph::Class()) {
1514  TGraph *g = (TGraph*)fPoly;
1515  by = g->GetY();
1516  bn = g->GetN();
1517  fYmin = by[0];
1518  for (i=1; i<bn; i++) {if (fYmin > by[i]) fYmin = by[i];}
1519  }
1520 
1521  if (fPoly->IsA() == TMultiGraph::Class()) {
1522  TMultiGraph *mg = (TMultiGraph*)fPoly;
1523  TList *gl = mg->GetListOfGraphs();
1524  if (!gl) return fYmin;
1525  TGraph *g;
1526  TIter next(gl);
1527  Bool_t first = kTRUE;
1528  while ((g = (TGraph*) next())) {
1529  by = g->GetY();
1530  bn = g->GetN();
1531  if (first) {fYmin = by[0]; first = kFALSE;}
1532  for (i=0; i<bn; i++) {if (fYmin > by[i]) fYmin = by[i];}
1533  }
1534  }
1535 
1536  return fYmin;
1537 }
1538 
1539 ////////////////////////////////////////////////////////////////////////////////
1540 /// Return "true" if the point (x,y) is inside the bin.
1541 
1542 Bool_t TH2PolyBin::IsInside(Double_t x, Double_t y) const
1543 {
1544  Int_t in=0;
1545 
1546  if (fPoly->IsA() == TGraph::Class()) {
1547  TGraph *g = (TGraph*)fPoly;
1548  in = g->IsInside(x, y);
1549  }
1550 
1551  if (fPoly->IsA() == TMultiGraph::Class()) {
1552  TMultiGraph *mg = (TMultiGraph*)fPoly;
1553  in = mg->IsInside(x, y);
1554  }
1555 
1556  return in;
1557 }
1558 
1559 ////////////////////////////////////////////////////////////////////////
1560 /// RE-implement dummy functions to avoid users calling the
1561 /// corresponding implemntations in TH1 or TH2
1562 //////////////////////////////////////////////////////////////////////////
1563 
1564 ////////////////////////////////////////////////////////////////////////////////
1565 /// Performs the operation: this = this + c1*f1. NOT IMPLEMENTED for TH2Poly
1566 Bool_t TH2Poly::Add(TF1 *, Double_t, Option_t *)
1567 {
1568  Error("Add","Not implement for TH2Poly");
1569  return kFALSE;
1570 }
1571 
1572 ////////////////////////////////////////////////////////////////////////////////
1573 /// Replace contents of this histogram by the addition of h1 and h2. NOT IMPLEMENTED for TH2Poly
1574 Bool_t TH2Poly::Add(const TH1 *, const TH1 *, Double_t, Double_t)
1575 {
1576  Error("Add","Not implement for TH2Poly");
1577  return kFALSE;
1578 }
1579 
1580 ////////////////////////////////////////////////////////////////////////////////
1581 /// Performs the operation: this = this / c1*f1. NOT IMPLEMENTED for TH2Poly
1582 Bool_t TH2Poly::Divide(TF1 *, Double_t)
1583 {
1584  Error("Divide","Not implement for TH2Poly");
1585  return kFALSE;
1586 }
1587 
1588 ////////////////////////////////////////////////////////////////////////////////
1589 /// NOT IMPLEMENTED for TH2Poly
1590 Bool_t TH2Poly::Multiply(TF1 *, Double_t)
1591 {
1592  Error("Multiply","Not implement for TH2Poly");
1593  return kFALSE;
1594 }
1595 ////////////////////////////////////////////////////////////////////////////////
1596 /// NOT IMPLEMENTED for TH2Poly
1597 Double_t TH2Poly::ComputeIntegral(Bool_t )
1598 {
1599  Error("ComputeIntegral","Not implement for TH2Poly");
1600  return TMath::QuietNaN();
1601 }
1602 ////////////////////////////////////////////////////////////////////////////////
1603 /// NOT IMPLEMENTED for TH2Poly
1604 TH1 * TH2Poly::FFT(TH1*, Option_t * )
1605 {
1606  Error("FFT","Not implement for TH2Poly");
1607  return nullptr;
1608 }
1609 ////////////////////////////////////////////////////////////////////////////////
1610 /// NOT IMPLEMENTED for TH2Poly
1611 TH1 * TH2Poly::GetAsymmetry(TH1* , Double_t, Double_t)
1612 {
1613  Error("GetAsymmetry","Not implement for TH2Poly");
1614  return nullptr;
1615 }
1616 ////////////////////////////////////////////////////////////////////////////////
1617 /// NOT IMPLEMENTED for TH2Poly
1618 Double_t TH2Poly::Interpolate(Double_t, Double_t)
1619 {
1620  Error("Interpolate","Not implement for TH2Poly");
1621  return TMath::QuietNaN();
1622 }