library: libHist #include "TProfile3D.h" |
TProfile3D
class description - header file - source file - inheritance tree (.pdf)
private:
Double_t* GetB()
Double_t* GetW()
Double_t* GetW2()
protected:
virtual Int_t BufferFill(Double_t, Double_t)
virtual Int_t BufferFill(Double_t, Double_t, Double_t)
virtual Int_t BufferFill(Double_t, Double_t, Double_t, Double_t)
virtual Int_t BufferFill(Double_t x, Double_t y, Double_t z, Double_t t, Double_t w)
public:
TProfile3D()
TProfile3D(const char* name, const char* title, Int_t nbinsx, Double_t xlow, Double_t xup, Int_t nbinsy, Double_t ylow, Double_t yup, Int_t nbinsz, Double_t zlow, Double_t zup, Option_t* option = "")
TProfile3D(const char* name, const char* title, Int_t nbinsx, const Double_t* xbins, Int_t nbinsy, const Double_t* ybins, Int_t nbinsz, const Double_t* zbins, Option_t* option = "")
TProfile3D(const TProfile3D& profile)
virtual ~TProfile3D()
virtual void Add(TF1* h1, Double_t c1 = 1, Option_t* option = "")
virtual void Add(const TH1* h1, Double_t c1 = 1)
virtual void Add(const TH1* h1, const TH1* h2, Double_t c1 = 1, Double_t c2 = 1)
static void Approximate(Bool_t approx = kTRUE)
virtual Int_t BufferEmpty(Int_t action = 0)
void BuildOptions(Double_t tmin, Double_t tmax, Option_t* option)
static TClass* Class()
virtual void Copy(TObject& hnew) const
virtual void Divide(TF1* h1, Double_t c1 = 1)
virtual void Divide(const TH1* h1)
virtual void Divide(const TH1* h1, const TH1* h2, Double_t c1 = 1, Double_t c2 = 1, Option_t* option = "")
virtual TH1* DrawCopy(Option_t* option = "") const
virtual Int_t Fill(Double_t)
virtual Int_t Fill(const char*, Double_t)
virtual Int_t Fill(Double_t, Double_t)
virtual Int_t Fill(Double_t, Double_t, Double_t)
virtual Int_t Fill(const char*, const char*, const char*, Double_t)
virtual Int_t Fill(const char*, Double_t, const char*, Double_t)
virtual Int_t Fill(const char*, const char*, Double_t, Double_t)
virtual Int_t Fill(Double_t, const char*, const char*, Double_t)
virtual Int_t Fill(Double_t, const char*, Double_t, Double_t)
virtual Int_t Fill(Double_t, Double_t, const char*, Double_t)
virtual Int_t Fill(Double_t x, Double_t y, Double_t z, Double_t t)
virtual Int_t Fill(Double_t x, Double_t y, Double_t z, Double_t t, Double_t w)
virtual Double_t GetBinContent(Int_t bin) const
virtual Double_t GetBinContent(Int_t, Int_t) const
virtual Double_t GetBinContent(Int_t binx, Int_t biny, Int_t binz) const
virtual Double_t GetBinEntries(Int_t bin) const
virtual Double_t GetBinError(Int_t bin) const
virtual Double_t GetBinError(Int_t, Int_t) const
virtual Double_t GetBinError(Int_t binx, Int_t biny, Int_t binz) const
Option_t* GetErrorOption() const
virtual void GetStats(Double_t* stats) const
virtual Double_t GetTmax() const
virtual Double_t GetTmin() const
virtual TClass* IsA() const
virtual Long64_t Merge(TCollection* list)
virtual void Multiply(TF1* h1, Double_t c1 = 1)
virtual void Multiply(const TH1* h1)
virtual void Multiply(const TH1* h1, const TH1* h2, Double_t c1 = 1, Double_t c2 = 1, Option_t* option = "")
TH3D* ProjectionXYZ(const char* name = "_pxyz", Option_t* option = "e") const
virtual void PutStats(Double_t* stats)
virtual void RebinAxis(Double_t x, const char* ax)
virtual void Reset(Option_t* option = "")
virtual void SavePrimitive(ostream& out, Option_t* option = "")
virtual void Scale(Double_t c1 = 1)
virtual void SetBinEntries(Int_t bin, Double_t w)
virtual void SetBins(Int_t, Double_t, Double_t)
virtual void SetBins(Int_t, const Double_t*)
virtual void SetBins(Int_t, Double_t, Double_t, Int_t, Double_t, Double_t)
virtual void SetBins(Int_t nbinsx, Double_t xmin, Double_t xmax, Int_t nbinsy, Double_t ymin, Double_t ymax, Int_t nbinsz, Double_t zmin, Double_t zmax)
virtual void SetBins(Int_t, const Double_t*, Int_t, const Double_t*)
virtual void SetBuffer(Int_t buffersize, Option_t* opt = "")
virtual void SetErrorOption(Option_t* option = "")
virtual void ShowMembers(TMemberInspector& insp, char* parent)
virtual void Streamer(TBuffer& b)
void StreamerNVirtual(TBuffer& b)
protected:
TArrayD fBinEntries number of entries per bin
EErrorType fErrorMode Option to compute errors
Double_t fTmin Lower limit in T (if set)
Double_t fTmax Upper limit in T (if set)
Bool_t fScaling !True when TProfile3D::Scale is called
Double_t fTsumwt Total Sum of weight*T
Double_t fTsumwt2 Total Sum of weight*T*T
static Bool_t fgApproximate bin error approximation option
______________________________________________________________________________
Profile3D histograms are used to display the mean
value of T and its RMS for each cell in X,Y,Z.
Profile3D histograms are in many cases an
The inter-relation of three measured quantities X, Y, Z and T can always
be visualized by a four-dimensional histogram or scatter-plot;
its representation on the line-printer is not particularly
satisfactory, except for sparse data. If T is an unknown (but single-valued)
approximate function of X,Y,Z this function is displayed by a profile3D histogram with
much better precision than by a scatter-plot.
The following formulae show the cumulated contents (capital letters) and the values
displayed by the printing or plotting routines (small letters) of the elements for cell I, J.
2
H(I,J,K) = sum T E(I,J,K) = sum T
l(I,J,K) = sum l L(I,J,K) = sum l
h(I,J,K) = H(I,J,K)/L(I,J,K) s(I,J,K) = sqrt(E(I,J,K)/L(I,J,K)- h(I,J,K)**2)
e(I,J,K) = s(I,J,K)/sqrt(L(I,J,K))
In the special case where s(I,J,K) is zero (eg, case of 1 entry only in one cell)
e(I,J,K) is computed from the average of the s(I,J,K) for all cells.
This simple/crude approximation was suggested in order to keep the cell
during a fit operation.
Example of a profile3D histogram
{
TCanvas *c1 = new TCanvas("c1","Profile histogram example",200,10,700,500);
hprof3d = new TProfile3D("hprof3d","Profile of pt versus px, py and pz",40,-4,4,40,-4,4,40,0,20);
Double_t px, py, pz, pt;
TRandom3 r(0);
for ( Int_t i=0; i<25000; i++) {
r.Rannor(px,py);
pz = px*px + py*py;
pt = r.Landau(0,1);
hprof3d->Fill(px,py,pz,pt,1);
}
hprof3d->Draw();
}
NOTE: A TProfile3D is drawn as it was a simple TH3
TProfile3D()
*-*-*-*-*-*Default constructor for Profile3D histograms*-*-*-*-*-*-*-*-*
*-* ============================================
~TProfile3D()
*-*-*-*-*-*Default destructor for Profile3D histograms*-*-*-*-*-*-*-*-*
*-* ===========================================
TProfile3D(const char *name,const char *title,Int_t nx,Double_t xlow,Double_t xup,Int_t ny,Double_t ylow,Double_t yup,Int_t nz, Double_t zlow,Double_t zup,Option_t *option)
*-*-*-*-*-*Normal Constructor for Profile histograms*-*-*-*-*-*-*-*-*-*
*-* ==========================================
The first eleven parameters are similar to TH3D::TH3D.
All values of t are accepted at filling time.
To fill a profile3D histogram, one must use TProfile3D::Fill function.
Note that when filling the profile histogram the function Fill
checks if the variable t is betyween fTmin and fTmax.
If a minimum or maximum value is set for the T scale before filling,
then all values below tmin or above tmax will be discarded.
Setting the minimum or maximum value for the T scale before filling
has the same effect as calling the special TProfile3D constructor below
where tmin and tmax are specified.
H(I,J,K) is printed as the cell contents. The errors computed are s(I,J,K) if CHOPT='S'
(spread option), or e(I,J,K) if CHOPT=' ' (error on mean).
See TProfile3D::BuildOptions for explanation of parameters
see other constructors below with all possible combinations of
fix and variable bin size like in TH3D.
void BuildOptions(Double_t tmin, Double_t tmax, Option_t *option)
*-*-*-*-*-*-*Set Profile3D histogram structure and options*-*-*-*-*-*-*-*-*
*-* =============================================
If a cell has N data points all with the same value T (especially
possible when dealing with integers), the spread in T for that cell
is zero, and the uncertainty assigned is also zero, and the cell is
ignored in making subsequent fits. If SQRT(T) was the correct error
in the case above, then SQRT(T)/SQRT(N) would be the correct error here.
In fact, any cell with non-zero number of entries N but with zero spread
should have an uncertainty SQRT(T)/SQRT(N).
Now, is SQRT(T)/SQRT(N) really the correct uncertainty?
that it is only in the case where the T variable is some sort
of counting statistics, following a Poisson distribution. This should
probably be set as the default case. However, T can be any variable
from an original NTUPLE, not necessarily distributed "Poissonly".
The computation of errors is based on the parameter option:
option:
' ' (Default) Errors are Spread/SQRT(N) for Spread.ne.0. ,
" " SQRT(T)/SQRT(N) for Spread.eq.0,N.gt.0 ,
" " 0. for N.eq.0
's' Errors are Spread for Spread.ne.0. ,
" " SQRT(T) for Spread.eq.0,N.gt.0 ,
" " 0. for N.eq.0
'i' Errors are Spread/SQRT(N) for Spread.ne.0. ,
" " 1./SQRT(12.*N) for Spread.eq.0,N.gt.0 ,
" " 0. for N.eq.0
The third case above corresponds to Integer T values for which the
uncertainty is +-0.5, with the assumption that the probability that T
takes any value between T-0.5 and T+0.5 is uniform (the same argument
goes for T uniformly distributed between T and T+1); this would be
useful if T is an ADC measurement, for example. Other, fancier options
would be possible, at the cost of adding one more parameter to the PROFILE2D
For example, if all T variables are distributed according to some
known Gaussian of standard deviation Sigma, then:
'G' Errors are Spread/SQRT(N) for Spread.ne.0. ,
" " Sigma/SQRT(N) for Spread.eq.0,N.gt.0 ,
" " 0. for N.eq.0
For example, this would be useful when all T's are experimental quantities
measured with the same instrument with precision Sigma.
void Add(const TH1 *h1, Double_t c1)
Performs the operation: this = this + c1*h1
void Add(const TH1 *h1, const TH1 *h2, Double_t c1, Double_t c2)
*-*-*-*-*Replace contents of this profile3D by the addition of h1 and h2*-*-*
*-* ===============================================================
this = c1*h1 + c2*h2
void Approximate(Bool_t approx)
static function
set the fgApproximate flag. When the flag is true, the function GetBinError
will approximate the bin error with the average profile error on all bins
in the following situation only
- the number of bins in the profile3D is less than 10404 (eg 100x100x100)
- the bin number of entries is small ( <5)
- the estimated bin error is extremely small compared to the bin content
(see TProfile3D::GetBinError)
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 z, Double_t t, 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
fBuffer[4] = z of first entry
fBuffer[5] = t of first entry
void Copy(TObject &obj)
*-*-*-*-*-*-*-*Copy a Profile3D histogram to a new profile2D histogram*-*-*-*
*-* =======================================================
void Divide(const TH1 *h1)
*-*-*-*-*-*-*-*-*-*-*Divide this profile2D by h1*-*-*-*-*-*-*-*-*-*-*-*-*
*-* ===========================
this = this/h1
void Divide(const TH1 *h1, const TH1 *h2, Double_t c1, Double_t c2, Option_t *option)
*-*-*-*-*Replace contents of this profile2D by the division of h1 by h2*-*-*
*-* ==============================================================
this = c1*h1/(c2*h2)
TH1 * DrawCopy(Option_t *option)
*-*-*-*-*-*-*-*Draw a copy of this profile3D histogram*-*-*-*-*-*-*-*-*-*-*
*-* =======================================
Double_t GetBinContent(Int_t bin)
*-*-*-*-*-*-*Return bin content of a Profile3D histogram*-*-*-*-*-*-*-*-*
*-* ===========================================
Double_t GetBinEntries(Int_t bin)
*-*-*-*-*-*-*Return bin entries of a Profile3D histogram*-*-*-*-*-*-*-*-*
*-* ===========================================
Double_t GetBinError(Int_t bin)
*-*-*-*-*-*-*Return bin error of a Profile3D histogram*-*-*-*-*-*-*-*-*
Computing errors: A moving field
=================================
The computation of errors for a TProfile3D has evolved with the versions
of ROOT. The difficulty is in computing errors for bins with low statistics.
- prior to version 3.10, we had no special treatment of low statistic bins.
As a result, these bins had huge errors. The reason is that the
expression eprim2 is very close to 0 (rounding problems) or 0.
- The algorithm is modified/protected for the case
when a TProfile3D is projected (ProjectionX). The previous algorithm
generated a N^2 problem when projecting a TProfile3D with a large number of
bins (eg 100000).
- in version 3.10/02, a new static function TProfile::Approximate
is introduced to enable or disable (default) the approximation.
(see also comments in TProfile::GetBinError)
Option_t * GetErrorOption()
*-*-*-*-*-*-*-*-*-*Return option to compute profile2D errors*-*-*-*-*-*-*-*
*-* =========================================
void GetStats(Double_t *stats)
fill the array stats from the contents of this profile
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
stats[9] = sumwxz
stats[10]= sumwyz
stats[11]= sumwt
stats[12]= sumwt2
If no axis-subrange is specified (via TAxis::SetRange), the array stats
is simply a copy of the statistics quantities computed at filling time.
If a sub-range is specified, the function recomputes these quantities
from the bin contents in the current axis range.
Long64_t Merge(TCollection *li)
Merge all histograms in the collection in 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.
void Multiply(const TH1 *)
*-*-*-*-*-*-*-*-*-*-*Multiply this profile2D by h1*-*-*-*-*-*-*-*-*-*-*-*
*-* =============================
this = this*h1
void Multiply(const TH1 *, const TH1 *, Double_t, Double_t, Option_t *)
*-*-*-*-*Replace contents of this profile2D by multiplication of h1 by h2*-*
*-* ================================================================
this = (c1*h1)*(c2*h2)
TH3D * ProjectionXYZ(const char *name, Option_t *option)
*-*-*-*-*Project this profile3D into a 3-D histogram along X,Y,Z*-*-*-*-*-*-*
*-* =====================================================
The projection is always of the type TH3D.
if option "E" is specified, the errors are computed. (default)
if option "B" is specified, the content of bin of the returned histogram
will be equal to the GetBinEntries(bin) of the profile,
if option "C=E" the bin contents of the projection are set to the
bin errors of the profile
void PutStats(Double_t *stats)
Replace current statistics with the values in array stats
void Reset(Option_t *option)
*-*-*-*-*-*-*-*-*-*Reset contents of a Profile3D histogram*-*-*-*-*-*-*-*
*-* =======================================
void RebinAxis(Double_t x, const char* ax)
Profile histogram is resized along ax such that x is in the axis range.
The new axis limits are recomputed by doubling iteratively
the current axis range until the specified value x is within the limits.
The algorithm makes a copy of the histogram, then loops on all bins
of the old histogram to fill the rebinned histogram.
Takes into account errors (Sumw2) if any.
The bit kCanRebin must be set before invoking this function.
Ex: h->SetBit(TH1::kCanRebin);
void Scale(Double_t c1)
*-*-*-*-*Multiply this profile2D by a constant c1*-*-*-*-*-*-*-*-*
*-* ========================================
this = c1*this
This function uses the services of TProfile3D::Add
void SetBinEntries(Int_t bin, Double_t w)
*-*-*-*-*-*-*-*-*Set the number of entries in bin*-*-*-*-*-*-*-*-*-*-*-*
*-* ================================
void SetErrorOption(Option_t *option)
*-*-*-*-*-*-*-*-*-*Set option to compute profile2D errors*-*-*-*-*-*-*-*
*-* =======================================
The computation of errors is based on the parameter option:
option:
' ' (Default) Errors are Spread/SQRT(N) for Spread.ne.0. ,
" " SQRT(T)/SQRT(N) for Spread.eq.0,N.gt.0 ,
" " 0. for N.eq.0
's' Errors are Spread for Spread.ne.0. ,
" " SQRT(T) for Spread.eq.0,N.gt.0 ,
" " 0. for N.eq.0
'i' Errors are Spread/SQRT(N) for Spread.ne.0. ,
" " 1./SQRT(12.*N) for Spread.eq.0,N.gt.0 ,
" " 0. for N.eq.0
See TProfile3D::BuildOptions for explanation of all options
Author: Rene Brun 17/05/2006
Last update: root/hist:$Name: $:$Id: TProfile3D.cxx,v 1.4 2006/07/03 16:10:46 brun 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.