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