TProfile


class description - source file - inheritance tree

class TProfile : public TH1D

    private:
virtual Int_t Fill(Axis_t) virtual void FillN(Int_t, Axis_t*, Double_t*, Int_t) Double_t* GetB() Double_t* GetW() Double_t* GetW2() virtual void SetBins(Int_t, Double_t, Double_t, Int_t, Double_t, Double_t) virtual void SetBins(Int_t, Double_t, Double_t, Int_t, Double_t, Double_t, Int_t, Double_t, Double_t) public:
TProfile TProfile() TProfile TProfile(const char* name, const char* title, Int_t nbinsx, Axis_t xlow, Axis_t xup, Option_t* option) TProfile TProfile(const char* name, const char* title, Int_t nbinsx, Axis_t xlow, Axis_t xup, Axis_t ylow, Axis_t yup, Option_t* option) TProfile TProfile(const char* name, const char* title, Int_t nbinsx, Float_t* xbins, Option_t* option) TProfile TProfile(const char* name, const char* title, Int_t nbinsx, Double_t* xbins, Option_t* option) TProfile TProfile(const TProfile& profile) virtual void ~TProfile() virtual void Add(TF1* h1, Double_t c1 = 1) virtual void Add(TH1* h1, Double_t c1 = 1) virtual void Add(TH1* h1, TH1* h2, Double_t c1 = 1, Double_t c2 = 1) void BuildOptions(Double_t ymin, Double_t ymax, Option_t* option) static TClass* Class() virtual void Copy(TObject& hnew) virtual void Divide(TF1* h1, Double_t c1 = 1) virtual void Divide(TH1* h1) virtual void Divide(TH1* h1, TH1* h2, Double_t c1 = 1, Double_t c2 = 1, Option_t* option) virtual TH1* DrawCopy(Option_t* option) virtual Int_t Fill(Axis_t x, Axis_t y) virtual Int_t Fill(Axis_t x, Axis_t y, Stat_t w) virtual void FillN(Int_t ntimes, Axis_t* x, Axis_t* y, Double_t* w, Int_t stride = 1) virtual Stat_t GetBinContent(Int_t bin) const virtual Stat_t GetBinEntries(Int_t bin) const virtual Stat_t GetBinError(Int_t bin) const Option_t* GetErrorOption() const virtual Double_t GetYmax() const virtual Double_t GetYmin() const virtual TClass* IsA() const virtual void Multiply(TF1* h1, Double_t c1 = 1) virtual void Multiply(TH1* h1) virtual void Multiply(TH1* h1, TH1* h2, Double_t c1 = 1, Double_t c2 = 1, Option_t* option) TH1D* ProjectionX(const char* name = _px, Option_t* option = e) virtual TH1* Rebin(Int_t ngroup = 2, const char* newname) virtual void Reset(Option_t* option) virtual void Scale(Double_t c1 = 1) virtual void SetBinEntries(Int_t bin, Stat_t w) virtual void SetBins(Int_t nbins, Double_t xmin, Double_t xmax) virtual void SetErrorOption(Option_t* option) virtual void ShowMembers(TMemberInspector& insp, char* parent) virtual void Streamer(TBuffer& b) void StreamerNVirtual(TBuffer& b)

Data Members

protected:
TArrayD fBinEntries number of entries per bin EErrorType fErrorMode Option to compute errors Double_t fYmin Lower limit in Y (if set) Double_t fYmax Upper limit in Y (if set)

Class Description

  Profile histograms are used to display the mean
  value of Y and its RMS for each bin in X. Profile histograms are in many cases an
  elegant replacement of two-dimensional histograms : the inter-relation of two
  measured quantities X and Y can always be visualized by a two-dimensional
  histogram or scatter-plot; its representation on the line-printer is not particularly
  satisfactory, except for sparse data. If Y is an unknown (but single-valued)
  approximate function of X, this function is displayed by a profile 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 bin J.

                                                    2
      H(J)  =  sum Y                  E(J)  =  sum Y
      l(J)  =  sum l                  L(J)  =  sum l
      h(J)  =  H(J)/L(J)              s(J)  =  sqrt(E(J)/L(J)- h(J)**2)
      e(J)  =  s(J)/sqrt(L(J))

  In the special case where s(J) is zero (eg, case of 1 entry only in one bin)
  e(J) is computed from the average of the s(J) for all bins.
  This simple/crude approximation was suggested in order to keep the bin
  during a fit operation.

           Example of a profile histogram with its graphics output
{
  TCanvas *c1 = new TCanvas("c1","Profile histogram example",200,10,700,500);
  hprof  = new TProfile("hprof","Profile of pz versus px",100,-4,4,0,20);
  Float_t px, py, pz;
  for ( Int_t i=0; i<25000; i++) {
     gRandom->Rannor(px,py);
     pz = px*px + py*py;
     hprof->Fill(px,pz,1);
  }
  hprof->Draw();
}

/*

*/



TProfile() : TH1D()
*-*-*-*-*-*Default constructor for Profile histograms*-*-*-*-*-*-*-*-*
*-*        ==========================================

~TProfile()
*-*-*-*-*-*Default destructor for Profile histograms*-*-*-*-*-*-*-*-*
*-*        =========================================

TProfile(const char *name,const char *title,Int_t nbins,Axis_t xlow,Axis_t xup,Option_t *option) : TH1D(name,title,nbins,xlow,xup)
*-*-*-*-*-*Normal Constructor for Profile histograms*-*-*-*-*-*-*-*-*-*
*-*        ==========================================

  The first five parameters are similar to TH1D::TH1D.
  All values of y are accepted at filling time.
  To fill a profile histogram, one must use TProfile::Fill function.

  Note that when filling the profile histogram the function Fill
  checks if the variable y is betyween fYmin and fYmax.
  If a minimum or maximum value is set for the Y scale before filling,
  then all values below ymin or above ymax will be discarded.
  Setting the minimum or maximum value for the Y scale before filling
  has the same effect as calling the special TProfile constructor below
  where ymin and ymax are specified.

  H(J) is printed as the channel contents. The errors displayed are s(J) if CHOPT='S'
  (spread option), or e(J) if CHOPT=' ' (error on mean).

        See TProfile::BuildOptions for explanation of parameters


TProfile(const char *name,const char *title,Int_t nbins,Float_t *xbins,Option_t *option) : TH1D(name,title,nbins,xbins)
*-*-*-*-*-*Constructor for Profile histograms with variable bin size*-*-*-*-*
*-*        =========================================================

        See TProfile::BuildOptions for more explanations on errors


TProfile(const char *name,const char *title,Int_t nbins,Double_t *xbins,Option_t *option) : TH1D(name,title,nbins,xbins)
*-*-*-*-*-*Constructor for Profile histograms with variable bin size*-*-*-*-*
*-*        =========================================================

        See TProfile::BuildOptions for more explanations on errors


TProfile(const char *name,const char *title,Int_t nbins,Axis_t xlow,Axis_t xup,Axis_t ylow,Axis_t yup,Option_t *option) : TH1D(name,title,nbins,xlow,xup)
*-*-*-*-*-*Constructor for Profile histograms with range in y*-*-*-*-*-*
*-*        ==================================================
  The first five parameters are similar to TH1D::TH1D.
  Only the values of Y between YMIN and YMAX will be considered at filling time.
  ymin and ymax will also be the maximum and minimum values
  on the y scale when drawing the profile.

        See TProfile::BuildOptions for more explanations on errors


void BuildOptions(Double_t ymin, Double_t ymax, Option_t *option)
*-*-*-*-*-*-*Set Profile histogram structure and options*-*-*-*-*-*-*-*-*
*-*          ===========================================

    If a bin has N data points all with the same value Y (especially
    possible when dealing with integers), the spread in Y for that bin
    is zero, and the uncertainty assigned is also zero, and the bin is
    ignored in making subsequent fits. If SQRT(Y) was the correct error
    in the case above, then SQRT(Y)/SQRT(N) would be the correct error here.
    In fact, any bin with non-zero number of entries N but with zero spread
    should have an uncertainty SQRT(Y)/SQRT(N).

    Now, is SQRT(Y)/SQRT(N) really the correct uncertainty?
    that it is only in the case where the Y variable is some sort
    of counting statistics, following a Poisson distribution. This should
    probably be set as the default case. However, Y 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(Y)/SQRT(N) for Spread.eq.0,N.gt.0 ,
                      "     "  0.  for N.eq.0
     's'            Errors are Spread  for Spread.ne.0. ,
                      "     "  SQRT(Y)  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 Y values for which the
    uncertainty is +-0.5, with the assumption that the probability that Y
    takes any value between Y-0.5 and Y+0.5 is uniform (the same argument
    goes for Y uniformly distributed between Y and Y+1); this would be
    useful if Y is an ADC measurement, for example. Other, fancier options
    would be possible, at the cost of adding one more parameter to the PROFILE
    command. For example, if all Y 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 Y's are experimental quantities
    measured with the same instrument with precision Sigma.



TProfile(const TProfile &profile)

void Add(TF1 *, Double_t )
 Performs the operation: this = this + c1*f1

void Add(TH1 *h1, Double_t c1)
 Performs the operation: this = this + c1*h1

void Add(TH1 *h1, TH1 *h2, Double_t c1, Double_t c2)
*-*-*-*-*Replace contents of this profile by the addition of h1 and h2*-*-*
*-*      =============================================================

   this = c1*h1 + c2*h2


void Copy(TObject &obj)
*-*-*-*-*-*-*-*Copy a Profile histogram to a new profile histogram*-*-*-*-*
*-*            ===================================================

void Divide(TF1 *, Double_t )
 Performs the operation: this = this/(c1*f1)

void Divide(TH1 *h1)
*-*-*-*-*-*-*-*-*-*-*Divide this profile by h1*-*-*-*-*-*-*-*-*-*-*-*-*
*-*                  =========================

   this = this/h1


void Divide(TH1 *h1, TH1 *h2, Double_t c1, Double_t c2, Option_t *option)
*-*-*-*-*Replace contents of this profile by the division of h1 by h2*-*-*
*-*      ============================================================

   this = c1*h1/(c2*h2)


TH1* DrawCopy(Option_t *option)
*-*-*-*-*-*-*-*Draw a copy of this profile histogram*-*-*-*-*-*-*-*-*-*-*-*
*-*            =====================================

Int_t Fill(Axis_t x, Axis_t y)
*-*-*-*-*-*-*-*-*-*-*Fill a Profile histogram (no weights)*-*-*-*-*-*-*-*
*-*                  =====================================

Int_t Fill(Axis_t x, Axis_t y, Stat_t w)
*-*-*-*-*-*-*-*-*-*-*Fill a Profile histogram with weights*-*-*-*-*-*-*-*
*-*                  =====================================

void FillN(Int_t ntimes, Axis_t *x, Axis_t *y, Stat_t *w, Int_t stride)
*-*-*-*-*-*-*-*-*-*-*Fill a Profile histogram with weights*-*-*-*-*-*-*-*
*-*                  =====================================

Stat_t GetBinContent(Int_t bin) const
*-*-*-*-*-*-*Return bin content of a Profile histogram*-*-*-*-*-*-*-*-*-*
*-*          =========================================

Stat_t GetBinEntries(Int_t bin) const
*-*-*-*-*-*-*Return bin entries of a Profile histogram*-*-*-*-*-*-*-*-*-*
*-*          =========================================

Stat_t GetBinError(Int_t bin) const
*-*-*-*-*-*-*Return bin error of a Profile histogram*-*-*-*-*-*-*-*-*-*
*-*          =======================================

Option_t* GetErrorOption() const
*-*-*-*-*-*-*-*-*-*Return option to compute profile errors*-*-*-*-*-*-*-*-*
*-*                =======================================

void Multiply(TF1 *, Double_t )
 Performs the operation: this = this*c1*f1

void Multiply(TH1 *)
*-*-*-*-*-*-*-*-*-*-*Multiply this profile by h1*-*-*-*-*-*-*-*-*-*-*-*-*
*-*                  =============================

   this = this*h1


void Multiply(TH1 *, TH1 *, Double_t, Double_t, Option_t *)
*-*-*-*-*Replace contents of this profile by multiplication of h1 by h2*-*
*-*      ================================================================

   this = (c1*h1)*(c2*h2)


TH1D* ProjectionX(const char *name, Option_t *option)
*-*-*-*-*Project this profile into a 1-D histogram along X*-*-*-*-*-*-*
*-*      =================================================

   The projection is always of the type TH1D.

   if option "E" is specified, the errors are computed. (default)



TH1* Rebin(Int_t ngroup, const char*newname)
*-*-*-*-*Rebin this profile grouping ngroup bins together*-*-*-*-*-*-*-*-*
*-*      ================================================
   if newname is not blank a new temporary profile hnew is created.
   else the current profile is modified (default)
   The parameter ngroup indicates how many bins of this have to me merged
   into one bin of hnew
   If the original profile has errors stored (via Sumw2), the resulting
   profile has new errors correctly calculated.

   examples: if hp is an existing TProfile histogram with 100 bins
     hp->Rebin();  //merges two bins in one in hp: previous contents of hp are lost
     hp->Rebin(5); //merges five bins in one in hp
     TProfile *hnew = hp->Rebin(5,"hnew"); // creates a new profile hnew
                                       //merging 5 bins of hp in one bin

   NOTE1: If ngroup is not an exact divider of the number of bins,
          the top limit of the rebinned profile is changed
          to the upper edge of the bin=newbins*ngroup and the corresponding
          bins are added to the overflow bin.
          Statistics will be recomputed from the new bin contents.

void Reset(Option_t *option)
*-*-*-*-*-*-*-*-*-*Reset contents of a Profile histogram*-*-*-*-*-*-*-*-*
*-*                =====================================

void Scale(Double_t c1)
*-*-*-*-*Multiply this profile by a constant c1*-*-*-*-*-*-*-*-*
*-*      ======================================

   this = c1*this

 This function uses the services of TProfile::Add


void SetBinEntries(Int_t bin, Stat_t w)
*-*-*-*-*-*-*-*-*Set the number of entries in bin*-*-*-*-*-*-*-*-*-*-*-*
*-*              ================================

void SetBins(Int_t nx, Double_t xmin, Double_t xmax)
*-*-*-*-*-*-*-*-*Redefine  x axis parameters*-*-*-*-*-*-*-*-*-*-*-*
*-*              ===========================

void SetErrorOption(Option_t *option)
*-*-*-*-*-*-*-*-*-*Set option to compute profile errors*-*-*-*-*-*-*-*-*
*-*                =====================================

    The computation of errors is based on the parameter option:
    option:
     ' '  (Default) Errors are Spread/SQRT(N) for Spread.ne.0. ,
                      "     "  SQRT(Y)/SQRT(N) for Spread.eq.0,N.gt.0 ,
                      "     "  0.  for N.eq.0
     's'            Errors are Spread  for Spread.ne.0. ,
                      "     "  SQRT(Y)  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 TProfile::BuildOptions for explanation of all options

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



Inline Functions


               void SetBins(Int_t, Double_t, Double_t, Int_t, Double_t, Double_t, Int_t, Double_t, Double_t)
          Double_t* GetB()
          Double_t* GetW()
          Double_t* GetW2()
              Int_t Fill(Axis_t x, Axis_t y, Stat_t w)
               void FillN(Int_t ntimes, Axis_t* x, Axis_t* y, Double_t* w, Int_t stride = 1)
           Double_t GetYmin() const
           Double_t GetYmax() const
               void SetBins(Int_t nbins, Double_t xmin, Double_t xmax)
            TClass* Class()
            TClass* IsA() const
               void ShowMembers(TMemberInspector& insp, char* parent)
               void StreamerNVirtual(TBuffer& b)


Author: Rene Brun 29/09/95
Last update: root/hist:$Name: $:$Id: TProfile.cxx,v 1.9 2000/12/18 14:54:49 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.