library: libHist
#include "TH2.h"

TH2


class description - source file - inheritance tree (.pdf)

class TH2 : public TH1

Inheritance Chart:
TObject
<-
TNamed
TAttLine
TAttFill
TAttMarker
<-
TH1
<-
TH2
<-
TH2C
TH2D
<-
TProfile2D
TH2F
TH2I
TH2S

    protected:
virtual Int_t BufferFill(Double_t, Double_t) virtual Int_t BufferFill(Double_t x, Double_t y, Double_t w) public:
TH2() TH2(const char* name, const char* title, Int_t nbinsx, Double_t xlow, Double_t xup, Int_t nbinsy, Double_t ylow, Double_t yup) TH2(const char* name, const char* title, Int_t nbinsx, const Double_t* xbins, Int_t nbinsy, Double_t ylow, Double_t yup) TH2(const char* name, const char* title, Int_t nbinsx, Double_t xlow, Double_t xup, Int_t nbinsy, const Double_t* ybins) TH2(const char* name, const char* title, Int_t nbinsx, const Double_t* xbins, Int_t nbinsy, const Double_t* ybins) TH2(const char* name, const char* title, Int_t nbinsx, const Float_t* xbins, Int_t nbinsy, const Float_t* ybins) TH2(const TH2&) virtual ~TH2() virtual Int_t BufferEmpty(Int_t action = 0) virtual Double_t Chi2Test(const TH1* h, Option_t* option, Int_t constraint = 0) const static TClass* Class() virtual void Copy(TObject& hnew) const virtual Int_t Fill(Double_t) virtual Int_t Fill(const char*, Double_t) virtual Int_t Fill(Double_t x, Double_t y) virtual Int_t Fill(Double_t x, Double_t y, Double_t w) virtual Int_t Fill(Double_t x, const char* namey, Double_t w) virtual Int_t Fill(const char* namex, Double_t y, Double_t w) virtual Int_t Fill(const char* namex, const char* namey, Double_t w) virtual void FillN(Int_t, const Double_t*, const Double_t*, Int_t) virtual void FillN(Int_t ntimes, const Double_t* x, const Double_t* y, const Double_t* w, Int_t stride = 1) virtual void FillRandom(const char* fname, Int_t ntimes = 5000) virtual void FillRandom(TH1* h, Int_t ntimes = 5000) virtual void FitSlicesX(TF1* f1 = 0, Int_t binmin = 1, Int_t binmax = 0, Int_t cut = 0, Option_t* option = "QNR") virtual void FitSlicesY(TF1* f1 = 0, Int_t binmin = 1, Int_t binmax = 0, Int_t cut = 0, Option_t* option = "QNR") virtual Double_t GetBinWithContent2(Double_t c, Int_t& binx, Int_t& biny, Int_t firstx = 0, Int_t lastx = 0, Int_t firsty = 0, Int_t lasty = 0, Double_t maxdiff = 0) const virtual Double_t GetCorrelationFactor(Int_t axis1 = 1, Int_t axis2 = 2) const virtual Double_t GetCovariance(Int_t axis1 = 1, Int_t axis2 = 2) const virtual void GetRandom2(Double_t& x, Double_t& y) virtual void GetStats(Double_t* stats) const virtual Double_t Integral(Option_t* option = "") const virtual Double_t Integral(Int_t, Int_t, Option_t* = "") const virtual Double_t Integral(Int_t binx1, Int_t binx2, Int_t biny1, Int_t biny2, Option_t* option = "") const virtual Double_t Integral(Int_t, Int_t, Int_t, Int_t, Int_t, Int_t, Option_t* = "") const virtual TClass* IsA() const virtual Double_t KolmogorovTest(const TH1* h2, Option_t* option = "") const virtual Long64_t Merge(TCollection* list) TH2& operator=(const TH2&) TProfile* ProfileX(const char* name = "_pfx", Int_t firstybin = -1, Int_t lastybin = -1, Option_t* option = "") const TProfile* ProfileY(const char* name = "_pfy", Int_t firstxbin = -1, Int_t lastxbin = -1, Option_t* option = "") const TH1D* ProjectionX(const char* name = "_px", Int_t firstybin = -1, Int_t lastybin = -1, Option_t* option = "") const TH1D* ProjectionY(const char* name = "_py", Int_t firstxbin = -1, Int_t lastxbin = -1, Option_t* option = "") const virtual void PutStats(Double_t* stats) virtual TH2* Rebin2D(Int_t nxgroup = 2, Int_t nygroup = 2, const char* newname = "") virtual TH2* RebinX(Int_t ngroup = 2, const char* newname = "") virtual TH2* RebinY(Int_t ngroup = 2, const char* newname = "") virtual void Reset(Option_t* option = "") virtual void ShowMembers(TMemberInspector& insp, char* parent) virtual void Streamer(TBuffer& b) void StreamerNVirtual(TBuffer& b)

Data Members


    protected:
Double_t fScalefactor Scale factor Double_t fTsumwy Total Sum of weight*Y Double_t fTsumwy2 Total Sum of weight*Y*Y Double_t fTsumwxy Total Sum of weight*X*Y

Class Description

 Service class for 2-Dim histogram classes

  TH2C a 2-D histogram with one byte per cell (char)
  TH2S a 2-D histogram with two bytes per cell (short integer)
  TH2I a 2-D histogram with four bytes per cell (32 bits integer)
  TH2F a 2-D histogram with four bytes per cell (float)
  TH2D a 2-D histogram with eight bytes per cell (double)


TH2()

TH2(const char *name,const char *title,Int_t nbinsx,Double_t xlow,Double_t xup ,Int_t nbinsy,Double_t ylow,Double_t yup) :TH1(name,title,nbinsx,xlow,xup)
 see comments in the TH1 base class constructors

TH2(const char *name,const char *title,Int_t nbinsx,const Double_t *xbins ,Int_t nbinsy,Double_t ylow,Double_t yup) :TH1(name,title,nbinsx,xbins)
 see comments in the TH1 base class constructors

TH2(const char *name,const char *title,Int_t nbinsx,Double_t xlow,Double_t xup ,Int_t nbinsy,const Double_t *ybins) :TH1(name,title,nbinsx,xlow,xup)
 see comments in the TH1 base class constructors

TH2(const char *name,const char *title,Int_t nbinsx,const Double_t *xbins ,Int_t nbinsy,const Double_t *ybins) :TH1(name,title,nbinsx,xbins)
 see comments in the TH1 base class constructors

TH2(const char *name,const char *title,Int_t nbinsx,const Float_t *xbins ,Int_t nbinsy,const Float_t *ybins) :TH1(name,title,nbinsx,xbins)
 see comments in the TH1 base class constructors

TH2(const TH2 &h) : TH1()
 Copy constructor.
 The list of functions is not copied. (Use Clone if needed)

~TH2()

Int_t BufferEmpty(Int_t action)
 Fill histogram with all entries in the buffer.
 action = -1 histogram is reset and refilled from the buffer (called by THistPainter::Paint)
 action =  0 histogram is filled from the buffer
 action =  1 histogram is filled and buffer is deleted
             The buffer is automatically deleted when the number of entries
             in the buffer is greater than the number of entries in the histogram

Int_t BufferFill(Double_t x, Double_t y, Double_t w)
 accumulate arguments in buffer. When buffer is full, empty the buffer
 fBuffer[0] = number of entries in buffer
 fBuffer[1] = w of first entry
 fBuffer[2] = x of first entry
 fBuffer[3] = y of first entry

Double_t Chi2Test(const TH1 *h, Option_t *option, Int_t constraint) const
The Chi2 (Pearson's) test for differences between h and this histogram.
a small value of prob indicates a significant difference between the distributions

if the data was collected in such a way that the number of entries
in the first histogram is necessarily equal to the number of entries
in the second, the parameter _constraint_ must be made 1. Default is 0.
any additional constraints on the data lower the number of degrees of freedom
(i.e. increase constraint to more positive values) in accordance with
their number

/options:
  "O" : overflows included
  "U" : underflows included
  by default underflows and overflows are not included

  "P"        : print information about number of degrees of freedom and the value of chi2
  "Chi2"     : the function returns the Chisquare instead of the probability
  "Chi2/ndf" : the function returns the Chi2/ndf
  if none of the options "Chi2" or "Chi2/ndf" is specified, the function returns
  the Pearson test, ie probability.

void Copy(TObject &obj) const

Int_t Fill(Double_t x,Double_t y)
*-*-*-*-*-*-*-*-*-*-*Increment cell defined by x,y by 1*-*-*-*-*-*-*-*-*-*
*-*                  ==================================
*-*
*-* if x or/and y is less than the low-edge of the corresponding axis first bin,
*-*   the Underflow cell is incremented.
*-* if x or/and y is greater than the upper edge of corresponding axis last bin,
*-*   the Overflow cell is incremented.
*-*
*-* If the storage of the sum of squares of weights has been triggered,
*-* via the function Sumw2, then the sum of the squares of weights is incremented
*-* by 1in the cell corresponding to x,y.
*-*
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*

Int_t Fill(Double_t x, Double_t y, Double_t w)
*-*-*-*-*-*-*-*-*-*-*Increment cell defined by x,y by a weight w*-*-*-*-*-*
*-*                  ===========================================
*-*
*-* if x or/and y is less than the low-edge of the corresponding axis first bin,
*-*   the Underflow cell is incremented.
*-* if x or/and y is greater than the upper edge of corresponding axis last bin,
*-*   the Overflow cell is incremented.
*-*
*-* If the storage of the sum of squares of weights has been triggered,
*-* via the function Sumw2, then the sum of the squares of weights is incremented
*-* by w^2 in the cell corresponding to x,y.
*-*
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*

Int_t Fill(const char *namex, const char *namey, Double_t w)
 Increment cell defined by namex,namey by a weight w

 if x or/and y is less than the low-edge of the corresponding axis first bin,
   the Underflow cell is incremented.
 if x or/and y is greater than the upper edge of corresponding axis last bin,
   the Overflow cell is incremented.

 If the storage of the sum of squares of weights has been triggered,
 via the function Sumw2, then the sum of the squares of weights is incremented
 by w^2 in the cell corresponding to x,y.


Int_t Fill(const char *namex, Double_t y, Double_t w)
 Increment cell defined by namex,y by a weight w

 if x or/and y is less than the low-edge of the corresponding axis first bin,
   the Underflow cell is incremented.
 if x or/and y is greater than the upper edge of corresponding axis last bin,
   the Overflow cell is incremented.

 If the storage of the sum of squares of weights has been triggered,
 via the function Sumw2, then the sum of the squares of weights is incremented
 by w^2 in the cell corresponding to x,y.


Int_t Fill(Double_t x, const char *namey, Double_t w)
 Increment cell defined by x,namey by a weight w

 if x or/and y is less than the low-edge of the corresponding axis first bin,
   the Underflow cell is incremented.
 if x or/and y is greater than the upper edge of corresponding axis last bin,
   the Overflow cell is incremented.

 If the storage of the sum of squares of weights has been triggered,
 via the function Sumw2, then the sum of the squares of weights is incremented
 by w^2 in the cell corresponding to x,y.


void FillN(Int_t ntimes, const Double_t *x, const Double_t *y, const Double_t *w, Int_t stride)
*-*-*-*-*-*-*Fill a 2-D histogram with an array of values and weights*-*-*-*
*-*          ========================================================
*-*
*-* ntimes:  number of entries in arrays x and w (array size must be ntimes*stride)
*-* x:       array of x values to be histogrammed
*-* y:       array of y values to be histogrammed
*-* w:       array of weights
*-* stride:  step size through arrays x, y and w
*-*
*-* If the storage of the sum of squares of weights has been triggered,
*-* via the function Sumw2, then the sum of the squares of weights is incremented
*-* by w[i]^2 in the cell corresponding to x[i],y[i].
*-* if w is NULL each entry is assumed a weight=1
*-*
*-* NB: function only valid for a TH2x object
*-*
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*

void FillRandom(const char *fname, Int_t ntimes)
*-*-*-*-*-*-*Fill histogram following distribution in function fname*-*-*-*
*-*          =======================================================
*-*
*-*   The distribution contained in the function fname (TF2) is integrated
*-*   over the channel contents.
*-*   It is normalized to 1.
*-*   Getting one random number implies:
*-*     - Generating a random number between 0 and 1 (say r1)
*-*     - Look in which bin in the normalized integral r1 corresponds to
*-*     - Fill histogram channel
*-*   ntimes random numbers are generated
*-*
*-*  One can also call TF2::GetRandom2 to get a random variate from a function.
*-*
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-**-*-*-*-*-*-*-*

void FillRandom(TH1 *h, Int_t ntimes)
*-*-*-*-*-*-*Fill histogram following distribution in histogram h*-*-*-*
*-*          ====================================================
*-*
*-*   The distribution contained in the histogram h (TH2) is integrated
*-*   over the channel contents.
*-*   It is normalized to 1.
*-*   Getting one random number implies:
*-*     - Generating a random number between 0 and 1 (say r1)
*-*     - Look in which bin in the normalized integral r1 corresponds to
*-*     - Fill histogram channel
*-*   ntimes random numbers are generated
*-*
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-**-*-*-*-*-*-*-*

void FitSlicesX(TF1 *f1, Int_t binmin, Int_t binmax, Int_t cut, Option_t *option)
 Project slices along X in case of a 2-D histogram, then fit each slice
 with function f1 and make a histogram for each fit parameter
 Only bins along Y between binmin and binmax are considered.
 if f1=0, a gaussian is assumed
 Before invoking this function, one can set a subrange to be fitted along X
 via f1->SetRange(xmin,xmax)
 The argument option (default="QNR") can be used to change the fit options.
     "Q"  means Quiet mode
     "N"  means do not show the result of the fit
     "R"  means fit the function in the specified function range
     "G2" merge 2 consecutive bins along X
     "G3" merge 3 consecutive bins along X
     "G4" merge 4 consecutive bins along X
     "G5" merge 5 consecutive bins along X

 Note that the generated histograms are added to the list of objects
 in the current directory. It is the user's responsability to delete
 these histograms.

  Example: Assume a 2-d histogram h2
   Root > h2->FitSlicesX(); produces 4 TH1D histograms
          with h2_0 containing parameter 0(Constant) for a Gaus fit
                    of each bin in Y projected along X
          with h2_1 containing parameter 1(Mean) for a gaus fit
          with h2_2 containing parameter 2(RMS)  for a gaus fit
          with h2_chi2 containing the chisquare/number of degrees of freedom for a gaus fit

   Root > h2->FitSlicesX(0,15,22,10);
          same as above, but only for bins 15 to 22 along Y
          and only for bins in Y for which the corresponding projection
          along X has more than cut bins filled.

  NOTE: To access the generated histograms in the current directory, do eg:
     TH1D *h2_1 = (TH1D*)gDirectory->Get("h2_1");

void FitSlicesY(TF1 *f1, Int_t binmin, Int_t binmax, Int_t cut, Option_t *option)
 Project slices along Y in case of a 2-D histogram, then fit each slice
 with function f1 and make a histogram for each fit parameter
 Only bins along X between binmin and binmax are considered.
 if f1=0, a gaussian is assumed
 Before invoking this function, one can set a subrange to be fitted along Y
 via f1->SetRange(ymin,ymax)
 The argument option (default="QNR") can be used to change the fit options.
     "Q"  means Quiet mode
     "N"  means do not show the result of the fit
     "R"  means fit the function in the specified function range
     "G2" merge 2 consecutive bins along Y
     "G3" merge 3 consecutive bins along Y
     "G4" merge 4 consecutive bins along Y
     "G5" merge 5 consecutive bins along Y

 Note that the generated histograms are added to the list of objects
 in the current directory. It is the user's responsability to delete
 these histograms.

  Example: Assume a 2-d histogram h2
   Root > h2->FitSlicesY(); produces 4 TH1D histograms
          with h2_0 containing parameter 0(Constant) for a Gaus fit
                    of each bin in X projected along Y
          with h2_1 containing parameter 1(Mean) for a gaus fit
          with h2_2 containing parameter 2(RMS)  for a gaus fit
          with h2_chi2 containing the chisquare/number of degrees of freedom for a gaus fit

   Root > h2->FitSlicesY(0,15,22,10);
          same as above, but only for bins 15 to 22 along X
          and only for bins in X for which the corresponding projection
          along Y has more than cut bins filled.

  NOTE: To access the generated histograms in the current directory, do eg:
     TH1D *h2_1 = (TH1D*)gDirectory->Get("h2_1");

 A complete example of this function is given in 
tutorial:fitslicesy.C
 with the following output:
/* */

Double_t GetBinWithContent2(Double_t c, Int_t &binx, Int_t &biny, Int_t firstx, Int_t lastx, Int_t firsty, Int_t lasty, Double_t maxdiff) const
 compute first cell (binx,biny) in the range [firstx,lastx](firsty,lasty] for which
 diff = abs(cell_content-c) <= maxdiff
 In case several cells in the specified range with diff=0 are found
 the first cell found is returned in binx,biny.
 In case several cells in the specified range satisfy diff <=maxdiff
 the cell with the smallest difference is returned in binx,biny.
 In all cases the function returns the smallest difference.

 NOTE1: if firstx <= 0, firstx is set to bin 1
        if (lastx < firstx then firstx is set to the number of bins in X
        ie if firstx=0 and lastx=0 (default) the search is on all bins in X.
        if firsty <= 0, firsty is set to bin 1
        if (lasty < firsty then firsty is set to the number of bins in Y
        ie if firsty=0 and lasty=0 (default) the search is on all bins in Y.
 NOTE2: if maxdiff=0 (default), the first cell with content=c is returned.

Double_t GetCorrelationFactor(Int_t axis1, Int_t axis2) const
*-*-*-*-*-*-*-*Return correlation factor between axis1 and axis2*-*-*-*-*
*-*            ====================================================

Double_t GetCovariance(Int_t axis1, Int_t axis2) const
*-*-*-*-*-*-*-*Return covariance between axis1 and axis2*-*-*-*-*
*-*            ====================================================

void GetRandom2(Double_t &x, Double_t &y)
 return 2 random numbers along axis x and y distributed according
 the cellcontents of a 2-dim histogram

void GetStats(Double_t *stats) const
 fill the array stats from the contents of this histogram
 The array stats must be correctly dimensionned in the calling program.
 stats[0] = sumw
 stats[1] = sumw2
 stats[2] = sumwx
 stats[3] = sumwx2
 stats[4] = sumwy
 stats[5] = sumwy2
 stats[6] = sumwxy

Double_t Integral(Option_t *option) const
Return integral of bin contents. Only bins in the bins range are considered.
 By default the integral is computed as the sum of bin contents in the range.
 if option "width" is specified, the integral is the sum of
 the bin contents multiplied by the bin width in x and in y.

Double_t Integral(Int_t binx1, Int_t binx2, Int_t biny1, Int_t biny2, Option_t *option) const
Return integral of bin contents in range [binx1,binx2],[biny1,biny2]
 for a 2-D histogram
 By default the integral is computed as the sum of bin contents in the range.
 if option "width" is specified, the integral is the sum of
 the bin contents multiplied by the bin width in x and in y.

Double_t KolmogorovTest(const TH1 *h2, Option_t *option) const
  Statistical test of compatibility in shape between
  THIS histogram and h2, using Kolmogorov test.
     Default: Ignore under- and overflow bins in comparison

     option is a character string to specify options
         "U" include Underflows in test
         "O" include Overflows
         "N" include comparison of normalizations
         "D" Put out a line of "Debug" printout

   The returned function value is the probability of test
       (much less than one means NOT compatible)

  Code adapted by Rene Brun from original HBOOK routine HDIFF

Long64_t Merge(TCollection *list)
Add all histograms in the collection to this histogram.
This function computes the min/max for the axes,
compute a new number of bins, if necessary,
add bin contents, errors and statistics.
If overflows are present and limits are different the function will fail.
The function returns the total number of entries in the result histogram
if the merge is successfull, -1 otherwise.

IMPORTANT remark. The 2 axis x and y may have different number
of bins and different limits, BUT the largest bin width must be
a multiple of the smallest bin width and the upper limit must also
be a multiple of the bin width.

TH2* RebinX(Int_t ngroup, const char *newname)
 Rebin only the X axis
 see Rebin2D

TH2* RebinY(Int_t ngroup, const char *newname)
 Rebin only the Y axis
 see Rebin2D

TH2* Rebin2D(Int_t nxgroup, Int_t nygroup, const char *newname)
   -*-*-*Rebin this histogram grouping nxgroup/nygroup bins along the xaxis/yaxis together*-*-*-*-
         =================================================================================
   if newname is not blank a new temporary histogram hnew is created.
   else the current histogram is modified (default)
   The parameter nxgroup/nygroup indicate how many bins along the xaxis/yaxis of this
   have to me merged into one bin of hnew
   If the original histogram has errors stored (via Sumw2), the resulting
   histograms has new errors correctly calculated.

   examples: if hpxpy is an existing TH2 histogram with 40 x 40 bins
     hpxpy->Rebin();  // merges two bins along the xaxis and yaxis in one in hpxpy
                      // Carefull: previous contents of hpxpy are lost
     hpxpy->RebinX(5); //merges five bins along the xaxis in one in hpxpy
     TH2 *hnew = hpxpy->RebinY(5,"hnew"); // creates a new histogram hnew
                                          // merging 5 bins of h1 along the yaxis in one bin

   NOTE : If nxgroup/nygroup is not an exact divider of the number of bins,
          along the xaxis/yaxis the top limit(s) of the rebinned histogram
          is changed to the upper edge of the xbin=newxbins*nxgroup resp.
          ybin=newybins*nygroup and the corresponding bins are added to
          the overflow bin.
          Statistics will be recomputed from the new bin contents.

TProfile* ProfileX(const char *name, Int_t firstybin, Int_t lastybin, Option_t *option) const
*-*-*-*-*Project a 2-D histogram into a profile histogram along X*-*-*-*-*-*
*-*      ========================================================

   The projection is made from the channels along the Y axis
   ranging from firstybin to lastybin included.
   By default, bins 1 to ny are included
   When all bins are included, the number of entries in the projection
   is set to the number of entries of the 2-D histogram, otherwise
   the number of entries is incremented by 1 for all non empty cells.

   if option "d" is specified, the profile is drawn in the current pad.

   Using a TCutG object, it is possible to select a sub-range of a 2-D histogram.
   One must create a graphical cut (mouse or C++) and specify the name
   of the cut between [] in the option.
   For example, with a TCutG named "cutg", one can call:
      myhist->ProfileX(" ",firstybin,lastybin,"[cutg]");
   To invert the cut, it is enough to put a "-" in front of its name:
      myhist->ProfileX(" ",firstybin,lastybin,"[-cutg]");
   It is possible to apply several cuts ("," means logical AND):
      myhist->ProfileX(" ",firstybin,lastybin,[cutg1,cutg2]");

   NOTE that if a TProfile named name exists in the current directory or pad,
   the histogram is reset and filled again with the current contents of the TH2.
   The X axis attributes of the TH2 are copied to the X axis of the profile.

TProfile* ProfileY(const char *name, Int_t firstxbin, Int_t lastxbin, Option_t *option) const
*-*-*-*-*Project a 2-D histogram into a profile histogram along Y*-*-*-*-*-*
*-*      ========================================================

   The projection is made from the channels along the X axis
   ranging from firstxbin to lastxbin included.
   By default, bins 1 to nx are included
   When all bins are included, the number of entries in the projection
   is set to the number of entries of the 2-D histogram, otherwise
   the number of entries is incremented by 1 for all non empty cells.

   if option "d" is specified, the profile is drawn in the current pad.

   Using a TCutG object, it is possible to select a sub-range of a 2-D histogram.
   One must create a graphical cut (mouse or C++) and specify the name
   of the cut between [] in the option.
   For example, with a TCutG named "cutg", one can call:
      myhist->ProfileY(" ",firstybin,lastybin,"[cutg]");
   To invert the cut, it is enough to put a "-" in front of its name:
      myhist->ProfileY(" ",firstybin,lastybin,"[-cutg]");
   It is possible to apply several cuts:
      myhist->ProfileY(" ",firstybin,lastybin,[cutg1,cutg2]");

   NOTE that if a TProfile named name exists in the current directory or pad,
   the histogram is reset and filled again with the current contents of the TH2.
   The Y axis attributes of the TH2 are copied to the X axis of the profile.

TH1D* ProjectionX(const char *name, Int_t firstybin, Int_t lastybin, Option_t *option) const
*-*-*-*-*Project a 2-D histogram into a 1-D histogram along X*-*-*-*-*-*-*
*-*      ====================================================

   The projection is always of the type TH1D.
   The projection is made from the channels along the Y axis
   ranging from firstybin to lastybin included.
   By default, bins 1 to ny are included
   When all bins are included, the number of entries in the projection
   is set to the number of entries of the 2-D histogram, otherwise
   the number of entries is incremented by 1 for all non empty cells.

   To make the projection in X of the underflow bin in Y, use firstybin=lastybin=0;
   To make the projection in X of the overflow bin in Y, use firstybin=lastybin=ny+1;

   if option "e" is specified, the errors are computed.
   if option "d" is specified, the projection is drawn in the current pad.

   Using a TCutG object, it is possible to select a sub-range of a 2-D histogram.
   One must create a graphical cut (mouse or C++) and specify the name
   of the cut between [] in the option.
   For example, with a TCutG named "cutg", one can call:
      myhist->ProjectionX(" ",firstybin,lastybin,"[cutg]");
   To invert the cut, it is enough to put a "-" in front of its name:
      myhist->ProjectionX(" ",firstybin,lastybin,"[-cutg]");
   It is possible to apply several cuts:
      myhist->ProjectionX(" ",firstybin,lastybin,[cutg1,cutg2]");

   NOTE that if a TH1D named name exists in the current directory or pad,
   the histogram is reset and filled again with the current contents of the TH2.
   The X axis attributes of the TH2 are copied to the X axis of the projection.

TH1D* ProjectionY(const char *name, Int_t firstxbin, Int_t lastxbin, Option_t *option) const
*-*-*-*-*Project a 2-D histogram into a 1-D histogram along Y*-*-*-*-*-*-*
*-*      ====================================================

   The projection is always of the type TH1D.
   The projection is made from the channels along the X axis
   ranging from firstxbin to lastxbin included.
   By default, bins 1 to nx are included
   When all bins are included, the number of entries in the projection
   is set to the number of entries of the 2-D histogram, otherwise
   the number of entries is incremented by 1 for all non empty cells.

   To make the projection in Y of the underflow bin in X, use firstxbin=lastxbin=0;
   To make the projection in Y of the overflow bin in X, use firstxbin=lastxbin=nx+1;

   if option "e" is specified, the errors are computed.
   if option "d" is specified, the projection is drawn in the current pad.

   Using a TCutG object, it is possible to select a sub-range of a 2-D histogram.
   One must create a graphical cut (mouse or C++) and specify the name
   of the cut between [] in the option.
   For example, with a TCutG named "cutg", one can call:
      myhist->ProjectionY(" ",firstxbin,lastxbin,"[cutg]");
   To invert the cut, it is enough to put a "-" in front of its name:
      myhist->ProjectionY(" ",firstxbin,lastxbin,"[-cutg]");
   It is possible to apply several cuts:
      myhist->ProjectionY(" ",firstxbin,lastxbin,[cutg1,cutg2]");

   NOTE that if a TH1D named name exists in the current directory or pad,
   the histogram is reset and filled again with the current contents of the TH2.
   The Y axis attributes of the TH2 are copied to the X axis of the projection.

void PutStats(Double_t *stats)
 Replace current statistics with the values in array stats

void Reset(Option_t *option)
*-*-*-*-*-*-*-*Reset this histogram: contents, errors, etc*-*-*-*-*-*-*-*
*-*            ===========================================

void Streamer(TBuffer &R__b)
 Stream an object of class TH2.



Inline Functions


              Int_t BufferFill(Double_t x, Double_t y, Double_t w)
              Int_t Fill(const char* namex, Double_t y, Double_t w)
              Int_t Fill(const char* namex, const char* namey, Double_t w)
               void FillN(Int_t ntimes, const Double_t* x, const Double_t* y, const Double_t* w, Int_t stride = 1)
           Double_t Integral(Int_t binx1, Int_t binx2, Int_t biny1, Int_t biny2, Option_t* option = "") const
           Double_t Integral(Int_t, Int_t, Int_t, Int_t, Int_t, Int_t, Option_t* = "") const
            TClass* Class()
            TClass* IsA() const
               void ShowMembers(TMemberInspector& insp, char* parent)
               void StreamerNVirtual(TBuffer& b)
               TH2& operator=(const TH2&)


Author: Rene Brun 26/12/94
Last update: root/hist:$Name: $:$Id: TH2.cxx,v 1.87 2006/02/03 21:55:38 pcanal Exp $
Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *


ROOT page - Class index - Class Hierarchy - Top of the page

This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.