TH3
class description - source file - inheritance tree
class TH3 : public TH1, public TAtt3D
public:
TH3 TH3()
TH3 TH3(const char* name, const char* title, Int_t nbinsx, Axis_t xlow, Axis_t xup, Int_t nbinsy, Axis_t ylow, Axis_t yup, Int_t nbinsz, Axis_t zlow, Axis_t zup)
TH3 TH3(const char* name, const char* title, Int_t nbinsx, Float_t* xbins, Int_t nbinsy, Float_t* ybins, Int_t nbinsz, Float_t* zbins)
TH3 TH3(const char* name, const char* title, Int_t nbinsx, Double_t* xbins, Int_t nbinsy, Double_t* ybins, Int_t nbinsz, Double_t* zbins)
TH3 TH3(TH3&)
virtual void ~TH3()
static TClass* Class()
virtual void Copy(TObject& hnew)
virtual Int_t Fill(Axis_t)
virtual Int_t Fill(Axis_t, Stat_t)
virtual Int_t Fill(Axis_t x, Axis_t y, Axis_t z)
virtual Int_t Fill(Axis_t x, Axis_t y, Axis_t z, Stat_t w)
virtual void FillRandom(const char* fname, Int_t ntimes = 5000)
virtual void FillRandom(TH1* h, Int_t ntimes = 5000)
virtual void FitSlicesZ(TF1* f1 = 0, Int_t binminx = 1, Int_t binmaxx = 0, Int_t binminy = 1, Int_t binmaxy = 0, Int_t cut = 0, Option_t* option = QNR)
virtual void GetRandom3(Axis_t& x, Axis_t& y, Axis_t& z)
virtual void GetStats(Stat_t* stats) const
virtual Stat_t Integral()
virtual Stat_t Integral(Int_t, Int_t)
virtual Stat_t Integral(Int_t, Int_t, Int_t, Int_t)
virtual Stat_t Integral(Int_t binx1, Int_t binx2, Int_t biny1, Int_t biny2, Int_t binz1, Int_t binz2)
virtual TClass* IsA() const
virtual Double_t KolmogorovTest(TH1* h2, Option_t* option)
TH1* Project3D(Option_t* option = x)
TH1D* ProjectionZ(const char* name = _pz, Int_t firstxbin = 0, Int_t lastxbin = 9999, Int_t firstybin = 0, Int_t lastybin = 9999, Option_t* option)
virtual void PutStats(Stat_t* stats)
virtual void Reset(Option_t* option)
virtual void ShowMembers(TMemberInspector& insp, char* parent)
virtual void Sizeof3D() const
virtual void Streamer(TBuffer& b)
void StreamerNVirtual(TBuffer& b)
See also
-
TH3C, TH3D, TH3F, TH3S
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
*-*
*-* The 3-D histogram classes derived from the 1-D histogram classes.
*-* all operations are supported (fill, fit).
*-* Drawing is currently restricted to one single option.
*-* A cloud of points is drawn. The number of points is proportional to
*-* cell content.
*-*
TH3C a 3-D histogram with one byte per cell (char)
TH3S a 3-D histogram with two bytes per cell (short integer)
TH3F a 3-D histogram with four bytes per cell (float)
TH3D a 3-D histogram with eight bytes per cell (double)
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
TH3()
TH3(const char *name,const char *title,Int_t nbinsx,Axis_t xlow,Axis_t xup
,Int_t nbinsy,Axis_t ylow,Axis_t yup
,Int_t nbinsz,Axis_t zlow,Axis_t zup)
:TH1(name,title,nbinsx,xlow,xup),
TAtt3D()
*-*-*-*-*-*-*-*-*Normal constructor for fix bin size 3-D histograms*-*-*-*-*
*-* ==================================================
TH3(const char *name,const char *title,Int_t nbinsx,Float_t *xbins
,Int_t nbinsy,Float_t *ybins
,Int_t nbinsz,Float_t *zbins)
:TH1(name,title,nbinsx,xbins),
TAtt3D()
*-*-*-*-*-*-*-*Normal constructor for variable bin size 3-D histograms*-*-*-*
*-* =======================================================
TH3(const char *name,const char *title,Int_t nbinsx,Double_t *xbins
,Int_t nbinsy,Double_t *ybins
,Int_t nbinsz,Double_t *zbins)
:TH1(name,title,nbinsx,xbins),
TAtt3D()
*-*-*-*-*-*-*-*Normal constructor for variable bin size 3-D histograms*-*-*-*
*-* =======================================================
~TH3()
void Copy(TObject &obj)
Int_t Fill(Axis_t x, Axis_t y, Axis_t z)
*-*-*-*-*-*-*-*-*-*-*Increment cell defined by x,y,z by 1 *-*-*-*-*
*-* ====================================
*-*
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
Int_t Fill(Axis_t x, Axis_t y, Axis_t z, Stat_t w)
*-*-*-*-*-*-*-*-*-*-*Increment cell defined by x,y,z by a weight 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^2 in the cell corresponding to x,y,z.
*-*
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
void FillRandom(const char *fname, Int_t ntimes)
*-*-*-*-*-*-*Fill histogram following distribution in function fname*-*-*-*
*-* =======================================================
*-*
*-* The distribution contained in the function fname (TF1) 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 TF1::GetRandom 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 (TH3) 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 FitSlicesZ(TF1 *f1, Int_t binminx, Int_t binmaxx, Int_t binminy, Int_t binmaxy, Int_t cut, Option_t *option)
Project slices along Z in case of a 3-D histogram, then fit each slice
with function f1 and make a 2-d histogram for each fit parameter
Only cells in the bin range [binminx,binmaxx] and [binminy,binmaxy] are considered.
if f1=0, a gaussian is assumed
Before invoking this function, one can set a subrange to be fitted along Z
via f1->SetRange(zmin,zmax)
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
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 3-d histogram h3
Root > h3->FitSlicesZ(); produces 4 TH2D histograms
with h3_0 containing parameter 0(Constant) for a Gaus fit
of each cell in X,Y projected along Z
with h3_1 containing parameter 1(Mean) for a gaus fit
with h3_2 containing parameter 2(RMS) for a gaus fit
with h3_chi2 containing the chisquare/number of degrees of freedom for a gaus fit
Root > h3->Fit(0,15,22,0,0,10);
same as above, but only for bins 15 to 22 along X
and only for cells in X,Y for which the corresponding projection
along Z has more than cut bins filled.
NOTE: To access the generated histograms in the current directory, do eg:
TH2D *h3_1 = (TH2D*)gDirectory->Get("h3_1");
void GetRandom3(Axis_t &x, Axis_t &y, Axis_t &z)
return 3 random numbers along axis x , y and z distributed according
the cellcontents of a 3-dim histogram
void GetStats(Stat_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
stats[7] = sumwz
stats[8] = sumwz2
Stat_t Integral()
Return integral of bin contents. Only bins in the bins range are considered.
Stat_t Integral(Int_t binx1, Int_t binx2, Int_t biny1, Int_t biny2, Int_t binz1, Int_t binz2)
Return integral of bin contents in range [binx1,binx2],[biny1,biny2],[binz1,binz2]
for a 3-D histogram
Double_t KolmogorovTest(TH1 *h2, Option_t *option)
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)
WARNING !!!! THIS FUNCTION NOT YET TESTED
I started from TH2::KolmogorovTest, but changes are probably required
when invoking KolmogorovProb to take into account the 3rd dimension
It would be nice if a mathematician could look into this.
Code adapted by Rene Brun from original HBOOK routine HDIFF
TH1D* ProjectionZ(const char *name, Int_t ixmin, Int_t ixmax, Int_t iymin, Int_t iymax, Option_t *option)
*-*-*-*-*Project a 3-D histogram into a 1-D histogram along Z*-*-*-*-*-*-*
*-* ====================================================
The projection is always of the type TH1D.
The projection is made from the cells along the X axis
ranging from ixmin to ixmax and iymin to iymax included.
if option "E" is specified, the errors are computed.
code from Paola Collins & Hans Dijkstra
TH1* Project3D(Option_t *option)
Project a 3-d histogram into 1 or 2-d histograms depending on the
option parameter
option may contain a combination of the characters x,y,z,e
option = "x" return the x projection into a TH1D histogram
option = "y" return the y projection into a TH1D histogram
option = "z" return the z projection into a TH1D histogram
option = "xy" return the x versus y projection into a TH2D histogram
option = "yx" return the y versus x projection into a TH2D histogram
option = "xz" return the x versus z projection into a TH2D histogram
option = "zx" return the z versus x projection into a TH2D histogram
option = "yz" return the y versus z projection into a TH2D histogram
option = "zy" return the z versus y projection into a TH2D histogram
If option contains the string "e", errors are computed
The projection is made for the selected bins only.
To select a bin range along an axis, use TAxis::SetRange, eg
h3.GetYAxis()->SetRange(23,56);
void PutStats(Stat_t *stats)
Replace current statistics with the values in array stats
void Reset(Option_t *option)
*-*-*-*-*-*-*-*Reset this histogram: contents, errors, etc*-*-*-*-*-*-*-*
*-* ===========================================
void Sizeof3D() const
*-*-*-*-*-*-*Return total size of this 3-D shape with its attributes*-*-*
*-* ==========================================================
void Streamer(TBuffer &R__b)
Stream an object of class TH3.
Inline Functions
Int_t Fill(Axis_t x, Axis_t y, Axis_t z)
Int_t Fill(Axis_t x, Axis_t y, Axis_t z, Stat_t w)
Stat_t Integral(Int_t, Int_t, Int_t, Int_t)
Stat_t Integral(Int_t binx1, Int_t binx2, Int_t biny1, Int_t biny2, Int_t binz1, Int_t binz2)
TClass* Class()
TClass* IsA() const
void ShowMembers(TMemberInspector& insp, char* parent)
void StreamerNVirtual(TBuffer& b)
TH3 TH3(TH3&)
Author: Rene Brun 27/10/95
Last update: root/hist:$Name: $:$Id: TH3.cxx,v 1.8 2000/12/13 15:13:51 brun Exp $
Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
ROOT page - Class index - 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.