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