TH2Poly.cxx

Go to the documentation of this file.
00001 // @(#)root/hist:$Id: TH2Poly.cxx 37434 2010-12-09 12:56:54Z couet $
00002 // TH2Poly v2.1
00003 // Author: Olivier Couet, Deniz Gunceler
00004 
00005 /*************************************************************************
00006  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
00007  * All rights reserved.                                                  *
00008  *                                                                       *
00009  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00010  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00011  *************************************************************************/
00012 
00013 #include "TROOT.h"
00014 #include "TClass.h"
00015 #include "TH2Poly.h"
00016 #include "TCutG.h"
00017 #include "TList.h"
00018 #include "TMath.h"
00019 #include "TMultiGraph.h"
00020 #include "TGraph.h"
00021 #include "TStyle.h"
00022 #include "TCanvas.h"
00023 #include <stdlib.h>
00024 #include <string.h>
00025 #include <stdio.h>
00026 #include <ctype.h>
00027 #include "Riostream.h"
00028 
00029 ClassImp(TH2Poly)
00030 
00031 //______________________________________________________________________________
00032 /* Begin_Html
00033 <center><h2>TH2Poly: 2D Histogram with Polygonal Bins</h2></center>
00034 
00035 <h3>Overview</h3>
00036 <tt>TH2Poly</tt> is a 2D Histogram class (TH2) allowing to define polygonal
00037 bins of arbitary shape.
00038 <p>
00039 Each bin in the <tt>TH2Poly</tt> histogram is a <tt>TH2PolyBin</tt> object.
00040 <tt>TH2PolyBin</tt> is a very simple class containing the vertices (stored
00041 as <tt>TGraph</tt>s or <tt>TMultiGraph</tt>s ) and contents of the polygonal
00042 bin as well as several related functions.
00043 <p>
00044 Essentially, a <tt>TH2Poly</tt> is a TList of <tt>TH2PolyBin</tt> objects
00045 with methods to manipulate them.
00046 <p>
00047 Bins are defined using one of the <tt>AddBin()</tt> methods. The bin definition
00048 should be done before filling.
00049 <p>
00050 The histogram can be filled with <tt>Fill(Double_t x, Double_t y, Double_t w)
00051 </tt>. <tt>w</tt> is the weight.
00052 If no weight is specified, it is assumed to be 1.
00053 <p>
00054 Not all histogram's area need to be binned. Filling an area without bins,
00055 will falls into the overflows. Adding a bin is not retroactive; it doesn't
00056 affect previous fillings. A <tt>Fill()</tt> call, that
00057 was previously ignored due to the lack of a bin at the specified location, is
00058 not reconsidered when that location is binned later.
00059 <p>
00060 If there are two overlapping bins, the first one in the list will be incremented
00061 by <tt>Fill()</tt>.
00062 <p>
00063 The histogram may automatically extends its limits if a bin outside the
00064 histogram limits is added. This is done when the default constructor (with no
00065 arguments) is used. It generates a histogram with no limits along the X and Y
00066 axis. Adding bins to it will extend it up to a proper size.
00067 <p>
00068 <tt>TH2Poly</tt> implements a partitioning algorithm to speed up bins' filling.
00069 The partitioning algorithm divides the histogram into regions called cells.
00070 The bins that each cell intersects are recorded in an array of <tt>TList</tt>s.
00071 When a coordinate in the histogram is to be filled; the method (quickly) finds
00072 which cell the coordinate belongs.  It then only loops over the bins
00073 intersecting that cell to find the bin the input coordinate corresponds to.
00074 The partitioning of the histogram is updated continuously as each bin is added.
00075 The default number of cells on each axis is 25. This number could be set to
00076 another value in the constructor or adjusted later by calling the
00077 <tt>ChangePartition(Int_t, Int_t)</tt> method. The partitioning algorithm is
00078 considerably faster than the brute force algorithm (i.e. checking if each bin
00079 contains the input coordinates), especially if the histogram is to be filled
00080 many times.
00081 <p>
00082 The following very simple macro shows how to build and fill a <tt>TH2Poly</tt>:
00083 <pre>
00084 {
00085    TH2Poly *h2p = new TH2Poly();
00086 
00087    Double_t x1[] = {0, 5, 5};
00088    Double_t y1[] = {0, 0, 5};
00089    Double_t x2[] = {0, -1, -1, 0};
00090    Double_t y2[] = {0, 0, -1, -1};
00091    Double_t x3[] = {4, 3, 0, 1, 2.4};
00092    Double_t y3[] = {4, 3.7, 1, 4.7, 3.5};
00093 
00094    h2p->AddBin(3, x1, y1);
00095    h2p->AddBin(3, x2, y2);
00096    h2p->AddBin(3, x3, y3);
00097 
00098    h2p->Fill(   3,    1, 3); // fill bin 1
00099    h2p->Fill(-0.5, -0.5, 7); // fill bin 2
00100    h2p->Fill(-0.7, -0.5, 1); // fill bin 2
00101    h2p->Fill(   1,    3, 5); // fill bin 3
00102 }
00103 </pre>
00104 
00105 More examples can bin found in <tt>$ROOTSYS/tutorials/hist/th2poly*.C</tt>
00106 
00107 <h3>Partitioning Algorithm</h3>
00108 The partitioning algorithm forms an essential part of the <tt>TH2Poly</tt>
00109 class. It is implemented to speed up the filling of bins.
00110 <p>
00111 With the brute force approach, the filling is done in the following way:  An
00112 iterator loops over all bins in the <tt>TH2Poly</tt> and invokes the
00113 method <tt>IsInside()</tt> for each of them.
00114 This method checks if the input location is in that bin. If the filling
00115 coordinate is inside, the bin is filled. Looping over all the bin is
00116 very slow.
00117 <p>
00118 The alternative is to divide the histogram into virtual rectangular regions
00119 called "cells". Each cell stores the pointers of the bins intersecting it.
00120 When a coordinate is to be filled, the method finds which cell the coordinate
00121 falls into. Since the cells are rectangular, this can be done very quickly.
00122 It then only loops over the bins associated with that cell.
00123 <p>
00124 The addition of bins to the appropriate cells is done when the bin is added
00125 to the histogram. To do this, <tt>AddBin()</tt> calls the
00126 <tt>AddBinToPartition()</tt> method.
00127 This method adds the input bin to the partitioning matrix.
00128 <p>
00129 The number of partition cells per axis can be specified in the constructor.
00130 If it is not specified, the default value of 25 along each axis will be
00131 assigned. This value was chosen because it is small enough to avoid slowing
00132 down AddBin(), while being large enough to enhance Fill() by a considerable
00133 amount. Regardless of how it is initialized at construction time, it can be
00134 changed later with the <tt>ChangePartition()</tt> method.
00135 <tt>ChangePartition()</tt> deletes the
00136 old partition matrix and generates a new one with the specified number of cells
00137 on each axis.
00138 <p>
00139 The optimum number of partition cells per axis changes with the number of
00140 times <tt>Fill()</tt> will be called.  Although partitioning greatly speeds up
00141 filling, it also adds a constant time delay into the code. When <tt>Fill()</tt>
00142 is to be called many times, it is more efficient to divide the histogram into
00143 a large number cells. However, if the histogram is to be filled only a few
00144 times, it is better to divide into a small number of cells.
00145 End_Html */
00146 
00147 
00148 //______________________________________________________________________________
00149 TH2Poly::TH2Poly()
00150 {
00151    // Default Constructor. No boundaries specified.
00152 
00153    Initialize(0., 0., 0., 0., 25, 25);
00154    SetName("NoName");
00155    SetTitle("NoTitle");
00156    SetFloat();
00157 }
00158 
00159 
00160 //______________________________________________________________________________
00161 TH2Poly::TH2Poly(const char *name,const char *title, Double_t xlow,Double_t xup
00162                                              , Double_t ylow,Double_t yup)
00163 {
00164    // Constructor with specified name and boundaries,
00165    // but no partition cell number.
00166 
00167    Initialize(xlow, xup, ylow, yup, 25, 25);
00168    SetName(name);
00169    SetTitle(title);
00170    SetFloat(kFALSE);
00171 }
00172 
00173 
00174 //______________________________________________________________________________
00175 TH2Poly::TH2Poly(const char *name,const char *title,
00176            Int_t nX, Double_t xlow, Double_t xup,
00177            Int_t nY, Double_t ylow, Double_t yup)
00178 {
00179    // Constructor with specified name and boundaries and partition cell number.
00180 
00181    Initialize(xlow, xup, ylow, yup, nX, nY);
00182    SetName(name);
00183    SetTitle(title);
00184    SetFloat(kFALSE);
00185 }
00186 
00187 
00188 //______________________________________________________________________________
00189 TH2Poly::~TH2Poly()
00190 {
00191    // Destructor.
00192 }
00193 
00194 
00195 //______________________________________________________________________________
00196 Int_t TH2Poly::AddBin(TObject *poly)
00197 {
00198    // Adds a new bin to the histogram. It can be any object having the method
00199    // IsInside(). It returns the bin number in the histogram. It returns 0 if
00200    // it failed to add. To allow the histogram limits to expand when a bin
00201    // outside the limits is added, call SetFloat() before adding the bin.
00202 
00203    if (!poly) return 0;
00204 
00205    if (fBins == 0) {fBins = new TList();}
00206 
00207    fNcells++;
00208    TH2PolyBin *bin = new TH2PolyBin(poly, fNcells);
00209 
00210    // If the bin lies outside histogram boundaries, then extends the boundaries.
00211    // Also changes the partition information accordingly
00212    Bool_t flag = kFALSE;
00213    if (fFloat) {
00214       if (fXaxis.GetXmin() > bin->GetXMin()) {
00215          fXaxis.Set(100, bin->GetXMin(), fXaxis.GetXmax());
00216          flag = kTRUE;
00217       }
00218       if (fXaxis.GetXmax() < bin->GetXMax()) {
00219          fXaxis.Set(100, fXaxis.GetXmin(), bin->GetXMax());
00220          flag = kTRUE;
00221       }
00222       if (fYaxis.GetXmin() > bin->GetYMin()) {
00223          fYaxis.Set(100, bin->GetYMin(), fYaxis.GetXmax());
00224          flag = kTRUE;
00225       }
00226       if (fYaxis.GetXmax() < bin->GetYMax()) {
00227          fYaxis.Set(100, fYaxis.GetXmin(), bin->GetYMax());
00228          flag = kTRUE;
00229       }
00230       if (flag) ChangePartition(fCellX, fCellY);
00231    } else {
00232       /*Implement polygon clipping code here*/
00233    }
00234 
00235    fBins->Add((TObject*) bin);
00236    SetNewBinAdded(kTRUE);
00237 
00238    // Adds the bin to the partition matrix
00239    AddBinToPartition(bin);
00240 
00241    return fNcells;
00242 }
00243 
00244 
00245 //______________________________________________________________________________
00246 Int_t TH2Poly::AddBin(Int_t n, const Double_t *x, const Double_t *y)
00247 {
00248    // Adds a new bin to the histogram. The number of vertices and their (x,y)
00249    // coordinates are required as input. It returns the bin number in the
00250    // histogram.
00251 
00252    TGraph *g = new TGraph(n, x, y);
00253    Int_t bin = AddBin(g);
00254    return bin;
00255 }
00256 
00257 
00258 //______________________________________________________________________________
00259 Int_t TH2Poly::AddBin(Double_t x1, Double_t y1, Double_t x2, Double_t  y2)
00260 {
00261    // Add a new bin to the histogram. The bin shape is a rectangle.
00262    // It returns the bin number of the bin in the histogram.
00263 
00264    Double_t x[] = {x1, x1, x2, x2};
00265    Double_t y[] = {y1, y2, y2, y1};
00266    TGraph *g = new TGraph(4, x, y);
00267    Int_t bin = AddBin(g);
00268    return bin;
00269 }
00270 
00271 
00272 //______________________________________________________________________________
00273 void TH2Poly::AddBinToPartition(TH2PolyBin *bin)
00274 {
00275    // Adds the input bin into the partition cell matrix. This method is called
00276    // in AddBin() and ChangePartition().
00277 
00278    // Cell Info
00279    Int_t nl, nr, mb, mt; // Max/min indices of the cells that contain the bin
00280    Double_t xclipl, xclipr, yclipb, yclipt; // x and y coordinates of a cell
00281    Double_t binXmax, binXmin, binYmax, binYmin; // The max/min bin coordinates
00282 
00283    binXmax = bin->GetXMax();
00284    binXmin = bin->GetXMin();
00285    binYmax = bin->GetYMax();
00286    binYmin = bin->GetYMin();
00287    nl = (Int_t)(floor((binXmin - fXaxis.GetXmin())/fStepX));
00288    nr = (Int_t)(floor((binXmax - fXaxis.GetXmin())/fStepX));
00289    mb = (Int_t)(floor((binYmin - fYaxis.GetXmin())/fStepY));
00290    mt = (Int_t)(floor((binYmax - fYaxis.GetXmin())/fStepY));
00291 
00292    // Make sure the array indices are correct.
00293    if (nr>=fCellX) nr = fCellX-1;
00294    if (mt>=fCellY) mt = fCellY-1;
00295    if (nl<0)       nl = 0;
00296    if (mb<0)       mb = 0;
00297 
00298    fNCells = fCellX*fCellY;
00299 
00300    // Loop over all cells
00301    for (int i = nl; i <= nr; i++) {
00302       xclipl = fXaxis.GetXmin() + i*fStepX;
00303       xclipr = xclipl + fStepX;
00304       for (int j = mb; j <= mt; j++) {
00305          yclipb = fYaxis.GetXmin() + j*fStepY;
00306          yclipt = yclipb + fStepY;
00307 
00308          // If the bin is completely inside the cell,
00309          // add that bin to the cell then return
00310          if ((binXmin >= xclipl) && (binXmax <= xclipr) &&
00311              (binYmax <= yclipt) && (binYmin >= yclipb)){
00312             fCells[i + j*fCellX].Add((TObject*) bin);
00313             fIsEmpty[i + j*fCellX] = kFALSE;  // Makes the cell non-empty
00314             return;
00315          }
00316 
00317          // If any of the sides of the cell intersect with any side of the bin, add that bin then continue
00318          if (IsIntersecting(bin, xclipl, xclipr, yclipb, yclipt)) {
00319             fCells[i + j*fCellX].Add((TObject*) bin);
00320             fIsEmpty[i + j*fCellX] = kFALSE;  // Makes the cell non-empty
00321             continue;
00322          }
00323          // If a corner of the cell is inside the bin and since there is no intersection, then that cell completely inside the bin.
00324          if((bin->IsInside(xclipl,yclipb)) || (bin->IsInside(xclipl,yclipt))){
00325             fCells[i + j*fCellX].Add((TObject*) bin);
00326             fIsEmpty[i + j*fCellX] = kFALSE;  // Makes the cell non-empty
00327             fCompletelyInside[i + fCellX*j] = kTRUE;
00328             continue;
00329          }
00330          if((bin->IsInside(xclipr,yclipb)) || (bin->IsInside(xclipr,yclipt))){
00331             fCells[i + j*fCellX].Add((TObject*) bin);
00332             fIsEmpty[i + j*fCellX] = kFALSE;  // Makes the cell non-empty
00333             fCompletelyInside[i + fCellX*j] = kTRUE;
00334             continue;
00335          }
00336       }
00337    }
00338 }
00339 
00340 
00341 //______________________________________________________________________________
00342 void TH2Poly::ChangePartition(Int_t n, Int_t m)
00343 {
00344    // Changes the number of partition cells in the histogram.
00345    // Deletes the old partition and constructs a new one.
00346 
00347    fCellX = n;                          // Set the number of cells
00348    fCellY = m;                          // Set the number of cells
00349 
00350    delete [] fCells;                    // Deletes the old partition
00351 
00352    fNCells = fCellX*fCellY;
00353    fCells  = new TList [fNCells];  // Sets an empty partition
00354 
00355    fStepX = (fXaxis.GetXmax() - fXaxis.GetXmin())/fCellX;      // Calculate cell width
00356    fStepY = (fYaxis.GetXmax() - fYaxis.GetXmin())/fCellY;      // Calculate cell height
00357 
00358    delete [] fIsEmpty;
00359    delete [] fCompletelyInside;
00360    fIsEmpty = new Bool_t [fNCells];
00361    fCompletelyInside = new Bool_t [fNCells];
00362 
00363    // Initializes the flags
00364    for (int i = 0; i<fNCells; i++) {
00365       fIsEmpty[i]          = kTRUE;
00366       fCompletelyInside[i] = kFALSE;
00367    }
00368 
00369    // TList iterator
00370    TIter    next(fBins);
00371    TObject  *obj;
00372 
00373    while((obj = next())){   // Loop over bins and add them to the partition
00374       AddBinToPartition((TH2PolyBin*) obj);
00375    }
00376 }
00377 
00378 
00379 //______________________________________________________________________________
00380 void TH2Poly::ClearBinContents()
00381 {
00382    // Clears the contents of all bins in the histogram.
00383 
00384    TIter next(fBins);
00385    TObject *obj;
00386    TH2PolyBin *bin;
00387 
00388    // Clears the bin contents
00389    while ((obj = next())) {
00390       bin = (TH2PolyBin*) obj;
00391       bin->ClearContent();
00392    }
00393 
00394    // Clears the statistics
00395    fTsumw   = 0;
00396    fTsumwx  = 0;
00397    fTsumwx2 = 0;
00398    fTsumwy  = 0;
00399    fTsumwy2 = 0;
00400    fEntries = 0;
00401 }
00402 
00403 
00404 //______________________________________________________________________________
00405 TH1 *TH2Poly::DrawCopy(Option_t *option) const
00406 {
00407    // Draw copy.
00408 
00409    TString opt = option;
00410    opt.ToLower();
00411    if (gPad && !opt.Contains("same")) gPad->Clear();
00412    TH2Poly *newth2 = (TH2Poly*)Clone();
00413    newth2->SetDirectory(0);
00414    newth2->SetBit(kCanDelete);
00415    newth2->AppendPad(option);
00416    return newth2;
00417 }
00418 
00419 
00420 //______________________________________________________________________________
00421 Int_t TH2Poly::FindBin(Double_t x, Double_t y, Double_t)
00422 {
00423    // Returns the bin number of the bin at the given coordinate. -1 to -9 are the
00424    // overflow and underflow bins.  overflow bin -5 is the unbinned areas in
00425    // the histogram (also called "the sea"). The third parameter can be left
00426    // blank.
00427    // The overflow/underflow bins are:
00428    //
00429    // -1 | -2 | -3
00430    // -------------
00431    // -4 | -5 | -6
00432    // -------------
00433    // -7 | -8 | -9
00434    //
00435    // where -5 means is the "sea" bin (i.e. unbinned areas)
00436 
00437 
00438    // Checks for overflow/underflow
00439    Int_t overflow = 0;
00440    if      (y > fYaxis.GetXmax()) overflow += -1;
00441    else if (y > fYaxis.GetXmin()) overflow += -4;
00442    else                           overflow += -7;
00443    if      (x > fXaxis.GetXmax()) overflow += -2;
00444    else if (x > fXaxis.GetXmin()) overflow += -1;
00445    if (overflow != -5) return overflow;  // If there is an overflow/undeflow, returns that value
00446 
00447    // Finds the cell (x,y) coordinates belong to
00448    Int_t n = (Int_t)(floor((x-fXaxis.GetXmin())/fStepX));
00449    Int_t m = (Int_t)(floor((y-fYaxis.GetXmin())/fStepY));
00450 
00451    // Make sure the array indices are correct.
00452    if (n>=fCellX) n = fCellX-1;
00453    if (m>=fCellY) m = fCellY-1;
00454    if (n<0)       n = 0;
00455    if (m<0)       m = 0;
00456 
00457    if (fIsEmpty[n+fCellX*m]) return -5;  // If the cell is empty, then the point must be on "the sea"
00458 
00459    TH2PolyBin *bin;
00460 
00461    TIter next(&fCells[n+fCellX*m]);
00462    TObject *obj;
00463 
00464    // Search for the bin in the cell
00465    while ((obj=next())) {
00466       bin  = (TH2PolyBin*)obj;
00467       if (bin->IsInside(x,y)) return bin->GetBinNumber();
00468    }
00469 
00470    // If the search has not returned a bin, the point must be on "the sea"
00471    return -5;
00472 }
00473 
00474 
00475 //______________________________________________________________________________
00476 Int_t TH2Poly::Fill(Double_t x, Double_t y)
00477 {
00478    // Increment the bin containing (x,y) by 1.
00479    // Uses the partitioning algorithm.
00480 
00481    return Fill(x, y, 1.0);
00482 }
00483 
00484 
00485 //______________________________________________________________________________
00486 Int_t TH2Poly::Fill(Double_t x, Double_t y, Double_t w)
00487 {
00488    // Increment the bin containing (x,y) by w.
00489    // Uses the partitioning algorithm.
00490 
00491    if (fNcells==0) return 0;  // If there are no bins in the histogram, does nothing
00492    Int_t overflow = 0;
00493    if      (y > fYaxis.GetXmax()) overflow += -1;
00494    else if (y > fYaxis.GetXmin()) overflow += -4;
00495    else                           overflow += -7;
00496    if      (x > fXaxis.GetXmax()) overflow += -2;
00497    else if(x > fXaxis.GetXmin())  overflow += -1;
00498    if (overflow != -5) {
00499       fOverflow[-overflow - 1]++;
00500       return 0;
00501    }
00502 
00503    // Finds the cell (x,y) coordinates belong to
00504    Int_t n = (Int_t)(floor((x-fXaxis.GetXmin())/fStepX));
00505    Int_t m = (Int_t)(floor((y-fYaxis.GetXmin())/fStepY));
00506 
00507    // Make sure the array indices are correct.
00508    if (n>=fCellX) n = fCellX-1;
00509    if (m>=fCellY) m = fCellY-1;
00510    if (n<0)       n = 0;
00511    if (m<0)       m = 0;
00512 
00513    if (fIsEmpty[n+fCellX*m]) return 0;  // If the cell is empty, then does not do anything
00514 
00515    TH2PolyBin *bin;
00516    Int_t bi;
00517 
00518    TIter next(&fCells[n+fCellX*m]);
00519    TObject *obj;
00520 
00521    while ((obj=next())) {
00522       bin  = (TH2PolyBin*)obj;
00523       bi = bin->GetBinNumber()-1;
00524       if (bin->IsInside(x,y)) {
00525          bin->Fill(w);
00526 
00527          // Statistics
00528          fTsumw   = fTsumw + w;       // Increments the total number of content
00529          fTsumwx  = fTsumwx + w*x;    // Increments the weighted sum of x coordinates
00530          fTsumwx2 = fTsumwx2 + w*x*x; // Increments the weighted sum of squares of x coordinates
00531          fTsumwy  = fTsumwy + w*y;    // Increments the weighted sum of y coordinates
00532          fTsumwy2 = fTsumwy2 + w*y*y; // Increments the weighted sum of squares of y coordinates
00533          if (fSumw2.fN) fSumw2.fArray[bi] += w*w;
00534          fEntries++;                  // Increments the total number of times Fill() was called
00535 
00536          SetBinContentChanged(kTRUE);
00537 
00538          return bin->GetBinNumber();
00539       }
00540    }
00541 
00542    fOverflow[4]++;  // Increments the "sea".
00543    return 0;
00544 }
00545 
00546 
00547 //______________________________________________________________________________
00548 Int_t TH2Poly::Fill(const char* name, Double_t w)
00549 {
00550    // Increment the bin named "name" by w.
00551 
00552    TString sname(name);
00553 
00554    TIter    next(fBins);
00555    TObject  *obj;
00556    TH2PolyBin *bin;
00557 
00558    while ((obj = next())) {
00559       bin = (TH2PolyBin*) obj;
00560       if (!sname.CompareTo(bin->GetPolygon()->GetName())) {
00561          bin->Fill(w);
00562          fEntries++;
00563          SetBinContentChanged(kTRUE);
00564          return bin->GetBinNumber();
00565       }
00566    }
00567 
00568    return 0;
00569 }
00570 
00571 
00572 //______________________________________________________________________________
00573 void TH2Poly::FillN(Int_t ntimes, const Double_t* x, const Double_t* y,
00574                                const Double_t* w, Int_t stride)
00575 {
00576    // Fills a 2-D histogram with an array of values and weights.
00577    //
00578    // ntimes:  number of entries in arrays x and w
00579    //          (array size must be ntimes*stride)
00580    // x:       array of x values to be histogrammed
00581    // y:       array of y values to be histogrammed
00582    // w:       array of weights
00583    // stride:  step size through arrays x, y and w
00584 
00585    for (int i = 0; i < ntimes; i += stride) {
00586       Fill(x[i], y[i], w[i]);
00587    }
00588 }
00589 
00590 
00591 //______________________________________________________________________________
00592 Double_t TH2Poly::Integral(Option_t* option) const
00593 {
00594    // Returns the integral of bin contents.
00595    // By default the integral is computed as the sum of bin contents.
00596    // If option "width" or "area" is specified, the integral is the sum of
00597    // the bin contents multiplied by the area of the bin.
00598 
00599    TString opt = option;
00600    opt.ToLower();
00601 
00602    if ((opt.Contains("width")) || (opt.Contains("area"))) {
00603       Double_t w;  // Weight of each bin.  Equals to the area of the bin if options "width" or "area" are specified.
00604       Double_t integral = 0.;
00605 
00606       TIter    next(fBins);
00607       TObject *obj;
00608       TH2PolyBin *bin;
00609       while ((obj=next())) {
00610          bin       = (TH2PolyBin*) obj;
00611          w         = bin->GetArea();
00612          integral += w*(bin->GetContent());
00613       }
00614 
00615       return integral;
00616    } else {
00617       return fTsumw;
00618    }
00619 }
00620 
00621 
00622 //______________________________________________________________________________
00623 Double_t TH2Poly::GetBinContent(Int_t bin) const
00624 {
00625    // Returns the content of the input bin
00626    // For the overflow/underflow/sea bins:
00627    //
00628    // -1 | -2 | -3
00629    // ---+----+----
00630    // -4 | -5 | -6
00631    // ---+----+----
00632    // -7 | -8 | -9
00633    //
00634    // where -5 is the "sea" bin (i.e. unbinned areas)
00635 
00636    if (bin>fNcells || bin==0) return 0;
00637    if (bin<0) return fOverflow[-bin - 1];
00638    return ((TH2PolyBin*) fBins->At(bin-1))->GetContent();
00639 }
00640 
00641 
00642 //______________________________________________________________________________
00643 Double_t TH2Poly::GetBinError(Int_t bin) const
00644 {
00645    // Returns the value of error associated to bin number bin.
00646    // If the sum of squares of weights has been defined (via Sumw2),
00647    // this function returns the sqrt(sum of w2).
00648    // otherwise it returns the sqrt(contents) for this bin.
00649 
00650    if (bin < 0) bin = 0;
00651    if (bin > (fNcells)) return 0;
00652    if (fBuffer) ((TH1*)this)->BufferEmpty();
00653    if (fSumw2.fN) {
00654       Double_t err2 = fSumw2.fArray[bin-1];
00655       return TMath::Sqrt(err2);
00656    }
00657    Double_t error2 = TMath::Abs(GetBinContent(bin));
00658    return TMath::Sqrt(error2);
00659 }
00660 
00661 
00662 //______________________________________________________________________________
00663 const char *TH2Poly::GetBinName(Int_t bin) const
00664 {
00665    // Returns the bin name.
00666 
00667    if (bin > (fNcells))  return "";
00668    if (bin < 0)          return "";
00669    return ((TH2PolyBin*) fBins->At(bin-1))->GetPolygon()->GetName();
00670 }
00671 
00672 
00673 //______________________________________________________________________________
00674 const char *TH2Poly::GetBinTitle(Int_t bin) const
00675 {
00676    // Returns the bin title.
00677 
00678    if (bin > (fNcells))  return "";
00679    if (bin < 0)          return "";
00680    return ((TH2PolyBin*) fBins->At(bin-1))->GetPolygon()->GetTitle();
00681 }
00682 
00683 
00684 //______________________________________________________________________________
00685 Double_t TH2Poly::GetMaximum() const
00686 {
00687    // Returns the maximum value of the histogram.
00688 
00689    if (fNcells==0) return 0;
00690 
00691    TH2PolyBin  *b;
00692 
00693    TIter next(fBins);
00694    TObject *obj;
00695    Double_t max,c;
00696 
00697    max = ((TH2PolyBin*) next())->GetContent();
00698 
00699    while ((obj=next())) {
00700       b = (TH2PolyBin*)obj;
00701       c = b->GetContent();
00702       if (c>max) max = c;
00703    }
00704    return max;
00705 }
00706 
00707 
00708 //______________________________________________________________________________
00709 Double_t TH2Poly::GetMaximum(Double_t maxval) const
00710 {
00711    // Returns the maximum value of the histogram that is less than maxval.
00712 
00713    if (fNcells==0) return 0;
00714 
00715    TH2PolyBin  *b;
00716 
00717    TIter next(fBins);
00718    TObject *obj;
00719    Double_t max,c;
00720 
00721    max = ((TH2PolyBin*) next())->GetContent();
00722 
00723    while ((obj=next())) {
00724       b = (TH2PolyBin*)obj;
00725       c = b->GetContent();
00726       if (c>max && c<maxval) max=c;
00727    }
00728    return max;
00729 }
00730 
00731 
00732 //______________________________________________________________________________
00733 Double_t TH2Poly::GetMinimum() const
00734 {
00735    // Returns the minimum value of the histogram.
00736 
00737    if (fNcells==0) return 0;
00738 
00739    TH2PolyBin  *b;
00740 
00741    TIter next(fBins);
00742    TObject *obj;
00743    Double_t min,c;
00744 
00745    min = ((TH2PolyBin*) next())->GetContent();
00746 
00747    while ((obj=next())) {
00748       b = (TH2PolyBin*)obj;
00749       c = b->GetContent();
00750       if (c<min) min=c;
00751    }
00752    return min;
00753 }
00754 
00755 
00756 //______________________________________________________________________________
00757 Double_t TH2Poly::GetMinimum(Double_t minval) const
00758 {
00759    // Returns the minimum value of the histogram that is greater than minval.
00760 
00761    if (fNcells==0) return 0;
00762 
00763    TH2PolyBin  *b;
00764 
00765    TIter next(fBins);
00766    TObject *obj;
00767    Double_t min,c;
00768 
00769    min = ((TH2PolyBin*) next())->GetContent();
00770 
00771    while ((obj=next())) {
00772       b = (TH2PolyBin*)obj;
00773       c = b->GetContent();
00774       if (c<min && c>minval) min=c;
00775    }
00776    return min;
00777 }
00778 
00779 
00780 //______________________________________________________________________________
00781 void TH2Poly::Honeycomb(Double_t xstart, Double_t ystart, Double_t a,
00782                      Int_t k, Int_t s)
00783 {
00784    // Bins the histogram using a honeycomb structure
00785 
00786    // Add the bins
00787    Double_t numberOfHexagonsInTheRow;
00788    Double_t x[6], y[6];
00789    Double_t xloop, yloop, xtemp;
00790    xloop = xstart; yloop = ystart + a/2.0;
00791    for (int sCounter = 0; sCounter < s; sCounter++) {
00792 
00793       xtemp = xloop; // Resets the temp variable
00794 
00795       // Determine the number of hexagons in that row
00796       if(sCounter%2 == 0){numberOfHexagonsInTheRow = k;}
00797       else{numberOfHexagonsInTheRow = k - 1;}
00798 
00799       for (int kCounter = 0; kCounter <  numberOfHexagonsInTheRow; kCounter++) {
00800 
00801          // Go around the hexagon
00802          x[0] = xtemp;
00803          y[0] = yloop;
00804          x[1] = x[0];
00805          y[1] = y[0] + a;
00806          x[2] = x[1] + a*TMath::Sqrt(3)/2.0;
00807          y[2] = y[1] + a/2.0;
00808          x[3] = x[2] + a*TMath::Sqrt(3)/2.0;
00809          y[3] = y[1];
00810          x[4] = x[3];
00811          y[4] = y[0];
00812          x[5] = x[2];
00813          y[5] = y[4] - a/2.0;
00814 
00815          this->AddBin(6, x, y);
00816 
00817          // Go right
00818          xtemp += a*TMath::Sqrt(3);
00819       }
00820 
00821       // Increment the starting position
00822       if (sCounter%2 == 0) xloop += a*TMath::Sqrt(3)/2.0;
00823       else                 xloop -= a*TMath::Sqrt(3)/2.0;
00824       yloop += 1.5*a;
00825    }
00826 }
00827 
00828 
00829 //______________________________________________________________________________
00830 void TH2Poly::Initialize(Double_t xlow, Double_t xup,
00831                       Double_t ylow, Double_t yup, Int_t n, Int_t m)
00832 {
00833    // Initializes the TH2Poly object.  This method is called by the constructor.
00834 
00835    Int_t i;
00836    fDimension = 2;  //The dimesion of the histogram
00837 
00838    fBins   = 0;
00839    fNcells = 0;
00840 
00841    // Sets the boundaries of the histogram
00842    fXaxis.Set(100, xlow, xup);
00843    fYaxis.Set(100, ylow, yup);
00844 
00845    for (i=0; i<9; i++) fOverflow[i] = 0.;
00846 
00847    // Statistics
00848    fEntries = 0;   // The total number of entries
00849    fTsumw   = 0.;  // Total amount of content in the histogram
00850    fTsumwx  = 0.;  // Weighted sum of x coordinates
00851    fTsumwx2 = 0.;  // Weighted sum of the squares of x coordinates
00852    fTsumwy2 = 0.;  // Weighted sum of the squares of y coordinates
00853    fTsumwy  = 0.;  // Weighted sum of y coordinates
00854 
00855    fCellX = n; // Set the number of cells to default
00856    fCellY = m; // Set the number of cells to default
00857 
00858    fNCells = fCellX*fCellY;
00859    fCells  = new TList [fNCells];  // Sets an empty partition
00860    fStepX  = (fXaxis.GetXmax() - fXaxis.GetXmin())/fCellX;     // Calculate cell width
00861    fStepY  = (fYaxis.GetXmax() - fYaxis.GetXmin())/fCellY;     // Calculate cell height
00862 
00863    fIsEmpty = new Bool_t [fNCells];                            // Declares the 'empty partition' flag
00864    fCompletelyInside = new Bool_t [fNCells];                   // Declares the 'cell is completely inside bin' flag
00865 
00866    for (i = 0; i<fNCells; i++) {   // Initializes the flags
00867       fIsEmpty[i] = kTRUE;
00868       fCompletelyInside[i] = kFALSE;
00869    }
00870 
00871    // 3D Painter flags
00872    SetNewBinAdded(kFALSE);
00873    SetBinContentChanged(kFALSE);
00874 }
00875 
00876 
00877 //______________________________________________________________________________
00878 Bool_t TH2Poly::IsIntersecting(TH2PolyBin *bin,
00879                                Double_t xclipl, Double_t xclipr,
00880                                Double_t yclipb, Double_t yclipt)
00881 {
00882    // Returns kTRUE if the input bin is intersecting with the
00883    // input rectangle (xclipl, xclipr, yclipb, yclipt)
00884 
00885    Int_t     gn;
00886    Double_t *gx;
00887    Double_t *gy;
00888    Bool_t inter = kFALSE;
00889    TObject *poly = bin->GetPolygon();
00890 
00891    if (poly->IsA() == TGraph::Class()) {
00892       TGraph *g = (TGraph*)poly;
00893       gx = g->GetX();
00894       gy = g->GetY();
00895       gn = g->GetN();
00896       inter = IsIntersectingPolygon(gn, gx, gy, xclipl, xclipr, yclipb, yclipt);
00897    }
00898 
00899    if (poly->IsA() == TMultiGraph::Class()) {
00900       TMultiGraph *mg = (TMultiGraph*)poly;
00901       TList *gl = mg->GetListOfGraphs();
00902       if (!gl) return inter;
00903       TGraph *g;
00904       TIter next(gl);
00905       while ((g = (TGraph*) next())) {
00906          gx = g->GetX();
00907          gy = g->GetY();
00908          gn = g->GetN();
00909          inter = IsIntersectingPolygon(gn, gx, gy, xclipl, xclipr, yclipb, yclipt);
00910          if (inter) return inter;
00911       }
00912    }
00913 
00914    return inter;
00915 }
00916 
00917 //______________________________________________________________________________
00918 Bool_t TH2Poly::IsIntersectingPolygon(Int_t bn, Double_t *x, Double_t *y,
00919                                       Double_t xclipl, Double_t xclipr,
00920                                       Double_t yclipb, Double_t yclipt)
00921 {
00922    // Returns kTRUE if the input polygon (bn, x, y) is intersecting with the
00923    // input rectangle (xclipl, xclipr, yclipb, yclipt)
00924 
00925    Bool_t p0R, p0L, p0T, p0B, p0xM, p0yM, p1R, p1L, p1T, p1B, p1xM, p1yM, p0In, p1In;
00926 
00927    for (int counter = 0; counter < (bn-1); counter++) {
00928       // If both are on the same side, return kFALSE
00929       p0L = x[counter]     <= xclipl; // Point 0 is on the left
00930       p1L = x[counter + 1] <= xclipl; // Point 1 is on the left
00931       if (p0L && p1L) continue;
00932       p0R = x[counter]     >= xclipr; // Point 0 is on the right
00933       p1R = x[counter + 1] >= xclipr; // Point 1 is on the right
00934       if (p0R && p1R) continue;
00935       p0T = y[counter]     >= yclipt; // Point 0 is at the top
00936       p1T = y[counter + 1] >= yclipt; // Point 1 is at the top
00937       if (p0T && p1T) continue;
00938       p0B = y[counter]     <= yclipb; // Point 0 is at the bottom
00939       p1B = y[counter + 1] <= yclipb; // Point 1 is at the bottom
00940       if (p0B && p1B) continue;
00941 
00942       // Checks to see if any are inside
00943       p0xM = !p0R && !p0L; // Point 0 is inside along x
00944       p0yM = !p0T && !p0B; // Point 1 is inside along x
00945       p1xM = !p1R && !p1L; // Point 0 is inside along y
00946       p1yM = !p1T && !p1B; // Point 1 is inside along y
00947       p0In = p0xM && p0yM; // Point 0 is inside
00948       p1In = p1xM && p1yM; // Point 1 is inside
00949       if (p0In) {
00950          if (p1In) continue;
00951          return kTRUE;
00952       } else {
00953          if (p1In) return kTRUE;
00954       }
00955 
00956       // We know by now that the points are not in the same side and not inside.
00957 
00958       // Checks to see if they are opposite
00959 
00960       if (p0xM && p1xM) return kTRUE;
00961       if (p0yM && p1yM) return kTRUE;
00962 
00963       // We now know that the points are in different x and y indices
00964 
00965       Double_t xcoord[3], ycoord[3];
00966       xcoord[0] = x[counter];
00967       xcoord[1] = x[counter + 1];
00968       ycoord[0] = y[counter];
00969       ycoord[1] = y[counter + 1];
00970 
00971       if (p0L) {
00972          if(p1T){
00973             xcoord[2] = xclipl;
00974             ycoord[2] = yclipb;
00975             if((TMath::IsInside(xclipl, yclipt, 3, xcoord, ycoord)) ||
00976                (TMath::IsInside(xclipr, yclipb, 3, xcoord, ycoord))) continue;
00977             else return kTRUE;
00978          } else if (p1B) {
00979             xcoord[2] = xclipl;
00980             ycoord[2] = yclipt;
00981             if((TMath::IsInside(xclipl, yclipb, 3, xcoord, ycoord)) ||
00982                (TMath::IsInside(xclipr, yclipt, 3, xcoord, ycoord))) continue;
00983             else return kTRUE;
00984          } else { // p1yM
00985             if (p0T) {
00986                xcoord[2] = xclipl;
00987                ycoord[2] = yclipb;
00988                if (TMath::IsInside(xclipr, yclipt, 3, xcoord, ycoord)) continue;
00989                else return kTRUE;
00990             }
00991             if (p0B) {
00992                xcoord[2] = xclipl;
00993                ycoord[2] = yclipt;
00994                if (TMath::IsInside(xclipr, yclipb, 3, xcoord, ycoord)) continue;
00995                else return kTRUE;
00996             }
00997          }
00998       } else if (p0R) {
00999          if (p1T) {
01000             xcoord[2] = xclipl;
01001             ycoord[2] = yclipb;
01002             if ((TMath::IsInside(xclipr, yclipb, 3, xcoord, ycoord)) ||
01003                 (TMath::IsInside(xclipl, yclipt, 3, xcoord, ycoord))) continue;
01004             else return kTRUE;
01005          } else if (p1B) {
01006             xcoord[2] = xclipl;
01007             ycoord[2] = yclipt;
01008             if ((TMath::IsInside(xclipl, yclipb, 3, xcoord, ycoord)) ||
01009                 (TMath::IsInside(xclipr, yclipt, 3, xcoord, ycoord))) continue;
01010             else return kTRUE;
01011          } else{ // p1yM
01012             if (p0T) {
01013                xcoord[2] = xclipr;
01014                ycoord[2] = yclipb;
01015                if (TMath::IsInside(xclipl, yclipt, 3, xcoord, ycoord)) continue;
01016                else return kTRUE;
01017             }
01018             if (p0B) {
01019                xcoord[2] = xclipr;
01020                ycoord[2] = yclipt;
01021                if (TMath::IsInside(xclipl, yclipb, 3, xcoord, ycoord)) continue;
01022                else return kTRUE;
01023             }
01024          }
01025       }
01026    }
01027    return kFALSE;
01028 }
01029 
01030 
01031 //______________________________________________________________________________
01032 void TH2Poly::SavePrimitive(ostream &out, Option_t *option)
01033 {
01034    // Save primitive as a C++ statement(s) on output stream out
01035 
01036    out <<"   "<<endl;
01037    out <<"   "<< ClassName() <<" *";
01038 
01039    //histogram pointer has by default the histogram name.
01040    //however, in case histogram has no directory, it is safer to add a incremental suffix
01041    static Int_t hcounter = 0;
01042    TString histName = GetName();
01043    if (!fDirectory && !histName.Contains("Graph")) {
01044       hcounter++;
01045       histName += "__";
01046       histName += hcounter;
01047    }
01048    const char *hname = histName.Data();
01049 
01050    //Construct the class initialization
01051    out << hname << " = new " << ClassName() << "(\"" << hname << "\", \""
01052        << GetTitle() << "\", " << fCellX << ", " << fXaxis.GetXmin() << ", " << fXaxis.GetXmax()
01053        << ", " << fCellY << ", " << fYaxis.GetXmin() << ", " << fYaxis.GetXmax() << ");" << endl;
01054 
01055    // Save Bins
01056    TIter       next(fBins);
01057    TObject    *obj;
01058    TH2PolyBin *th2pBin;
01059 
01060    while((obj = next())){
01061       th2pBin = (TH2PolyBin*) obj;
01062       th2pBin->GetPolygon()->SavePrimitive(out, Form("th2poly%s",histName.Data()));
01063    }
01064 
01065    // save bin contents
01066    out<<"   "<<endl;
01067    Int_t bin;
01068    for (bin=1;bin<=fNcells;bin++) {
01069       Double_t bc = GetBinContent(bin);
01070       if (bc) {
01071          out<<"   "<<hname<<"->SetBinContent("<<bin<<","<<bc<<");"<<endl;
01072       }
01073    }
01074 
01075    // save bin errors
01076    if (fSumw2.fN) {
01077       for (bin=1;bin<=fNcells;bin++) {
01078          Double_t be = GetBinError(bin);
01079          if (be) {
01080             out<<"   "<<hname<<"->SetBinError("<<bin<<","<<be<<");"<<endl;
01081          }
01082       }
01083    }
01084    TH1::SavePrimitiveHelp(out, hname, option);
01085 }
01086 
01087 
01088 //______________________________________________________________________________
01089 void TH2Poly::SetBinContent(Int_t bin, Double_t content)
01090 {
01091    // Sets the contents of the input bin to the input content
01092 
01093    if (bin > (fNcells)) return;
01094    ((TH2PolyBin*) fBins->At(bin-1))->SetContent(content);
01095    SetBinContentChanged(kTRUE);
01096 }
01097 
01098 
01099 //______________________________________________________________________________
01100 void TH2Poly::SetFloat(Bool_t flag)
01101 {
01102    // When set to kTRUE, allows the histogram to expand if a bin outside the
01103    // limits is added.
01104 
01105    fFloat = flag;
01106 }
01107 
01108 
01109 //______________________________________________________________________________
01110 TH2PolyBin::TH2PolyBin()
01111 {
01112    // Default constructor.
01113 
01114    fPoly    = 0;
01115    fContent = 0.;
01116    fNumber  = 0;
01117    fXmax    = -1111;
01118    fXmin    = -1111;
01119    fYmax    = -1111;
01120    fYmin    = -1111;
01121    fArea    = 0;
01122    SetChanged(kTRUE);
01123 }
01124 
01125 
01126 //______________________________________________________________________________
01127 TH2PolyBin::TH2PolyBin(TObject *poly, Int_t bin_number)
01128 {
01129    // Normal constructor.
01130 
01131    fContent = 0.;
01132    fNumber  = bin_number;
01133    fArea    = 0.;
01134    fPoly    = poly;
01135    fXmax    = -1111;
01136    fXmin    = -1111;
01137    fYmax    = -1111;
01138    fYmin    = -1111;
01139    SetChanged(kTRUE);
01140 }
01141 
01142 
01143 //______________________________________________________________________________
01144 TH2PolyBin::~TH2PolyBin()
01145 {
01146    // Destructor.
01147 
01148    if (fPoly) delete fPoly;
01149 }
01150 
01151 
01152 //______________________________________________________________________________
01153 Double_t TH2PolyBin::GetArea()
01154 {
01155    // Returns the area of the bin.
01156 
01157    Int_t     bn;
01158 
01159    if (fArea == 0) {
01160       if (fPoly->IsA() == TGraph::Class()) {
01161          TGraph *g = (TGraph*)fPoly;
01162          bn    = g->GetN();
01163          fArea = g->Integral(0,bn-1);
01164       }
01165 
01166       if (fPoly->IsA() == TMultiGraph::Class()) {
01167          TMultiGraph *mg = (TMultiGraph*)fPoly;
01168          TList *gl = mg->GetListOfGraphs();
01169          if (!gl) return fArea;
01170          TGraph *g;
01171          TIter next(gl);
01172          while ((g = (TGraph*) next())) {
01173             bn    = g->GetN();
01174             fArea = fArea + g->Integral(0,bn-1);
01175          }
01176       }
01177    }
01178 
01179    return fArea;
01180 }
01181 
01182 
01183 //______________________________________________________________________________
01184 Double_t TH2PolyBin::GetXMax()
01185 {
01186    // Returns the maximum value for the x coordinates of the bin.
01187 
01188    if (fXmax != -1111) return fXmax;
01189 
01190    Int_t     bn,i;
01191    Double_t *bx;
01192 
01193    if (fPoly->IsA() == TGraph::Class()) {
01194       TGraph *g = (TGraph*)fPoly;
01195       bx    = g->GetX();
01196       bn    = g->GetN();
01197       fXmax = bx[0];
01198       for (i=1; i<bn; i++) {if (fXmax < bx[i]) fXmax = bx[i];}
01199    }
01200 
01201    if (fPoly->IsA() == TMultiGraph::Class()) {
01202       TMultiGraph *mg = (TMultiGraph*)fPoly;
01203       TList *gl = mg->GetListOfGraphs();
01204       if (!gl) return fXmax;
01205       TGraph *g;
01206       TIter next(gl);
01207       Bool_t first = kTRUE;
01208       while ((g = (TGraph*) next())) {
01209          bx = g->GetX();
01210          bn = g->GetN();
01211          if (first) {fXmax = bx[0]; first = kFALSE;}
01212          for (i=0; i<bn; i++) {if (fXmax < bx[i]) fXmax = bx[i];}
01213       }
01214    }
01215 
01216    return fXmax;
01217 }
01218 
01219 
01220 //______________________________________________________________________________
01221 Double_t TH2PolyBin::GetXMin()
01222 {
01223    // Returns the minimum value for the x coordinates of the bin.
01224 
01225    if (fXmin != -1111) return fXmin;
01226 
01227    Int_t     bn,i;
01228    Double_t *bx;
01229 
01230    if (fPoly->IsA() == TGraph::Class()) {
01231       TGraph *g = (TGraph*)fPoly;
01232       bx    = g->GetX();
01233       bn    = g->GetN();
01234       fXmin = bx[0];
01235       for (i=1; i<bn; i++) {if (fXmin > bx[i]) fXmin = bx[i];}
01236    }
01237 
01238    if (fPoly->IsA() == TMultiGraph::Class()) {
01239       TMultiGraph *mg = (TMultiGraph*)fPoly;
01240       TList *gl = mg->GetListOfGraphs();
01241       if (!gl) return fXmin;
01242       TGraph *g;
01243       TIter next(gl);
01244       Bool_t first = kTRUE;
01245       while ((g = (TGraph*) next())) {
01246          bx = g->GetX();
01247          bn = g->GetN();
01248          if (first) {fXmin = bx[0]; first = kFALSE;}
01249          for (i=0; i<bn; i++) {if (fXmin > bx[i]) fXmin = bx[i];}
01250       }
01251    }
01252 
01253    return fXmin;
01254 }
01255 
01256 
01257 //______________________________________________________________________________
01258 Double_t TH2PolyBin::GetYMax()
01259 {
01260    // Returns the maximum value for the y coordinates of the bin.
01261 
01262    if (fYmax != -1111) return fYmax;
01263 
01264    Int_t     bn,i;
01265    Double_t *by;
01266 
01267    if (fPoly->IsA() == TGraph::Class()) {
01268       TGraph *g = (TGraph*)fPoly;
01269       by    = g->GetY();
01270       bn    = g->GetN();
01271       fYmax = by[0];
01272       for (i=1; i<bn; i++) {if (fYmax < by[i]) fYmax = by[i];}
01273    }
01274 
01275    if (fPoly->IsA() == TMultiGraph::Class()) {
01276       TMultiGraph *mg = (TMultiGraph*)fPoly;
01277       TList *gl = mg->GetListOfGraphs();
01278       if (!gl) return fYmax;
01279       TGraph *g;
01280       TIter next(gl);
01281       Bool_t first = kTRUE;
01282       while ((g = (TGraph*) next())) {
01283          by = g->GetY();
01284          bn = g->GetN();
01285          if (first) {fYmax = by[0]; first = kFALSE;}
01286          for (i=0; i<bn; i++) {if (fYmax < by[i]) fYmax = by[i];}
01287       }
01288    }
01289 
01290    return fYmax;
01291 }
01292 
01293 
01294 //______________________________________________________________________________
01295 Double_t TH2PolyBin::GetYMin()
01296 {
01297    // Returns the minimum value for the y coordinates of the bin.
01298 
01299    if (fYmin != -1111) return fYmin;
01300 
01301    Int_t     bn,i;
01302    Double_t *by;
01303 
01304    if (fPoly->IsA() == TGraph::Class()) {
01305       TGraph *g = (TGraph*)fPoly;
01306       by    = g->GetY();
01307       bn    = g->GetN();
01308       fYmin = by[0];
01309       for (i=1; i<bn; i++) {if (fYmin > by[i]) fYmin = by[i];}
01310    }
01311 
01312    if (fPoly->IsA() == TMultiGraph::Class()) {
01313       TMultiGraph *mg = (TMultiGraph*)fPoly;
01314       TList *gl = mg->GetListOfGraphs();
01315       if (!gl) return fYmin;
01316       TGraph *g;
01317       TIter next(gl);
01318       Bool_t first = kTRUE;
01319       while ((g = (TGraph*) next())) {
01320          by = g->GetY();
01321          bn = g->GetN();
01322          if (first) {fYmin = by[0]; first = kFALSE;}
01323          for (i=0; i<bn; i++) {if (fYmin > by[i]) fYmin = by[i];}
01324       }
01325    }
01326 
01327    return fYmin;
01328 }
01329 
01330 
01331 //______________________________________________________________________________
01332 Bool_t TH2PolyBin::IsInside(Double_t x, Double_t y) const
01333 {
01334    // Return "true" if the point (x,y) is inside the bin.
01335 
01336    Int_t in=0;
01337 
01338    if (fPoly->IsA() == TGraph::Class()) {
01339       TGraph *g = (TGraph*)fPoly;
01340       in = g->IsInside(x, y);
01341    }
01342 
01343    if (fPoly->IsA() == TMultiGraph::Class()) {
01344       TMultiGraph *mg = (TMultiGraph*)fPoly;
01345       in = mg->IsInside(x, y);
01346    }
01347 
01348    return in;
01349 }

Generated on Tue Jul 5 14:24:20 2011 for ROOT_528-00b_version by  doxygen 1.5.1