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