TProfile.cxx

Go to the documentation of this file.
00001 // @(#)root/hist:$Id: TProfile.cxx 38022 2011-02-09 16:12:32Z moneta $
00002 // Author: Rene Brun   29/09/95
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 #include "TProfile.h"
00013 #include "TMath.h"
00014 #include "TF1.h"
00015 #include "THLimitsFinder.h"
00016 #include "Riostream.h"
00017 #include "TVirtualPad.h"
00018 #include "TError.h"
00019 #include "TClass.h"
00020 
00021 #include "TProfileHelper.h"
00022 
00023 Bool_t TProfile::fgApproximate = kFALSE;
00024 
00025 ClassImp(TProfile)
00026 
00027 //______________________________________________________________________________
00028 //
00029 //  Profile histograms are used to display the mean
00030 //  value of Y and its RMS for each bin in X. Profile histograms are in many cases an
00031 //  elegant replacement of two-dimensional histograms : the inter-relation of two
00032 //  measured quantities X and Y can always be visualized by a two-dimensional
00033 //  histogram or scatter-plot; its representation on the line-printer is not particularly
00034 //  satisfactory, except for sparse data. If Y is an unknown (but single-valued)
00035 //  approximate function of X, this function is displayed by a profile histogram with
00036 //  much better precision than by a scatter-plot.
00037 //
00038 //  The following formulae show the cumulated contents (capital letters) and the values
00039 //  displayed by the printing or plotting routines (small letters) of the elements for bin J.
00040 //
00041 //                                                    2
00042 //      H(J)  =  sum Y                  E(J)  =  sum Y
00043 //      l(J)  =  sum l                  L(J)  =  sum l
00044 //      h(J)  =  H(J)/L(J)              s(J)  =  sqrt(E(J)/L(J)- h(J)**2)
00045 //      e(J)  =  s(J)/sqrt(L(J))
00046 //
00047 //  In the special case where s(J) is zero (eg, case of 1 entry only in one bin)
00048 //  e(J) is computed from the average of the s(J) for all bins if the static function
00049 //  TProfile::Approximate has been called.
00050 //  This simple/crude approximation was suggested in order to keep the bin
00051 //  during a fit operation. But note that this approximation is not the default behaviour.
00052 //
00053 //           Example of a profile histogram with its graphics output
00054 //{
00055 //  TCanvas *c1 = new TCanvas("c1","Profile histogram example",200,10,700,500);
00056 //  hprof  = new TProfile("hprof","Profile of pz versus px",100,-4,4,0,20);
00057 //  Float_t px, py, pz;
00058 //  for ( Int_t i=0; i<25000; i++) {
00059 //     gRandom->Rannor(px,py);
00060 //     pz = px*px + py*py;
00061 //     hprof->Fill(px,pz,1);
00062 //  }
00063 //  hprof->Draw();
00064 //}
00065 //Begin_Html
00066 /*
00067 <img src="gif/profile.gif">
00068 */
00069 //End_Html
00070 //
00071 
00072 //______________________________________________________________________________
00073 TProfile::TProfile() : TH1D()
00074 {
00075 //*-*-*-*-*-*Default constructor for Profile histograms*-*-*-*-*-*-*-*-*
00076 //*-*        ==========================================
00077 
00078    BuildOptions(0,0,"");
00079 }
00080 
00081 //______________________________________________________________________________
00082 TProfile::~TProfile()
00083 {
00084 //*-*-*-*-*-*Default destructor for Profile histograms*-*-*-*-*-*-*-*-*
00085 //*-*        =========================================
00086 
00087 }
00088 
00089 //______________________________________________________________________________
00090 TProfile::TProfile(const char *name,const char *title,Int_t nbins,Double_t xlow,Double_t xup,Option_t *option)
00091     : TH1D(name,title,nbins,xlow,xup)
00092 {
00093 //*-*-*-*-*-*Normal Constructor for Profile histograms*-*-*-*-*-*-*-*-*-*
00094 //*-*        ==========================================
00095 //
00096 //  The first five parameters are similar to TH1D::TH1D.
00097 //  All values of y are accepted at filling time.
00098 //  To fill a profile histogram, one must use TProfile::Fill function.
00099 //
00100 //  Note that when filling the profile histogram the function Fill
00101 //  checks if the variable y is betyween fYmin and fYmax.
00102 //  If a minimum or maximum value is set for the Y scale before filling,
00103 //  then all values below ymin or above ymax will be discarded.
00104 //  Setting the minimum or maximum value for the Y scale before filling
00105 //  has the same effect as calling the special TProfile constructor below
00106 //  where ymin and ymax are specified.
00107 //
00108 //  H(J) is printed as the channel contents. The errors displayed are s(J) if CHOPT='S'
00109 //  (spread option), or e(J) if CHOPT=' ' (error on mean).
00110 //
00111 //        See TProfile::BuildOptions for explanation of parameters
00112 //
00113 // see also comments in the TH1 base class constructors
00114 
00115    BuildOptions(0,0,option);
00116 }
00117 
00118 //______________________________________________________________________________
00119 TProfile::TProfile(const char *name,const char *title,Int_t nbins,const Float_t *xbins,Option_t *option)
00120     : TH1D(name,title,nbins,xbins)
00121 {
00122 //*-*-*-*-*-*Constructor for Profile histograms with variable bin size*-*-*-*-*
00123 //*-*        =========================================================
00124 //
00125 //        See TProfile::BuildOptions for more explanations on errors
00126 //
00127 // see also comments in the TH1 base class constructors
00128 
00129    BuildOptions(0,0,option);
00130 }
00131 
00132 //______________________________________________________________________________
00133 TProfile::TProfile(const char *name,const char *title,Int_t nbins,const Double_t *xbins,Option_t *option)
00134     : TH1D(name,title,nbins,xbins)
00135 {
00136 //*-*-*-*-*-*Constructor for Profile histograms with variable bin size*-*-*-*-*
00137 //*-*        =========================================================
00138 //
00139 //        See TProfile::BuildOptions for more explanations on errors
00140 //
00141 // see also comments in the TH1 base class constructors
00142 
00143    BuildOptions(0,0,option);
00144 }
00145 
00146 //______________________________________________________________________________
00147 TProfile::TProfile(const char *name,const char *title,Int_t nbins,const Double_t *xbins,Double_t ylow,Double_t yup,Option_t *option)
00148     : TH1D(name,title,nbins,xbins)
00149 {
00150 //*-*-*-*-*-*Constructor for Profile histograms with variable bin size*-*-*-*-*
00151 //*-*        =========================================================
00152 //
00153 //        See TProfile::BuildOptions for more explanations on errors
00154 //
00155 // see also comments in the TH1 base class constructors
00156 
00157    BuildOptions(ylow,yup,option);
00158 }
00159 
00160 //______________________________________________________________________________
00161 TProfile::TProfile(const char *name,const char *title,Int_t nbins,Double_t xlow,Double_t xup,Double_t ylow,Double_t yup,Option_t *option)
00162     : TH1D(name,title,nbins,xlow,xup)
00163 {
00164 //*-*-*-*-*-*Constructor for Profile histograms with range in y*-*-*-*-*-*
00165 //*-*        ==================================================
00166 //  The first five parameters are similar to TH1D::TH1D.
00167 //  Only the values of Y between ylow and yup will be considered at filling time.
00168 //  ylow and yup will also be the maximum and minimum values
00169 //  on the y scale when drawing the profile.
00170 //
00171 //        See TProfile::BuildOptions for more explanations on errors
00172 //
00173 // see also comments in the TH1 base class constructors
00174 
00175    BuildOptions(ylow,yup,option);
00176 }
00177 
00178 
00179 //______________________________________________________________________________
00180 void TProfile::BuildOptions(Double_t ymin, Double_t ymax, Option_t *option)
00181 {
00182 //*-*-*-*-*-*-*Set Profile histogram structure and options*-*-*-*-*-*-*-*-*
00183 //*-*          ===========================================
00184 //
00185 //    If a bin has N data points all with the same value Y (especially
00186 //    possible when dealing with integers), the spread in Y for that bin
00187 //    is zero, and the uncertainty assigned is also zero, and the bin is
00188 //    ignored in making subsequent fits. If SQRT(Y) was the correct error
00189 //    in the case above, then SQRT(Y)/SQRT(N) would be the correct error here.
00190 //    In fact, any bin with non-zero number of entries N but with zero spread
00191 //    should have an uncertainty SQRT(Y)/SQRT(N).
00192 //
00193 //    Now, is SQRT(Y)/SQRT(N) really the correct uncertainty?
00194 //    that it is only in the case where the Y variable is some sort
00195 //    of counting statistics, following a Poisson distribution. This should
00196 //    probably be set as the default case. However, Y can be any variable
00197 //    from an original NTUPLE, not necessarily distributed "Poissonly".
00198 //    The computation of errors is based on the parameter option:
00199 //    option:
00200 //     ' '  (Default) Errors are Spread/SQRT(N) for Spread.ne.0. ,
00201 //                      "     "  SQRT(Y)/SQRT(N) for Spread.eq.0,N.gt.0 ,
00202 //                      "     "  0.  for N.eq.0
00203 //     's'            Errors are Spread  for Spread.ne.0. ,
00204 //                      "     "  SQRT(Y)  for Spread.eq.0,N.gt.0 ,
00205 //                      "     "  0.  for N.eq.0
00206 //     'i'            Errors are Spread/SQRT(N) for Spread.ne.0. ,
00207 //                      "     "  1./SQRT(12.*N) for Spread.eq.0,N.gt.0 ,
00208 //                      "     "  0.  for N.eq.0
00209 //
00210 //    The third case above corresponds to Integer Y values for which the
00211 //    uncertainty is +-0.5, with the assumption that the probability that Y
00212 //    takes any value between Y-0.5 and Y+0.5 is uniform (the same argument
00213 //    goes for Y uniformly distributed between Y and Y+1); this would be
00214 //    useful if Y is an ADC measurement, for example. 
00215 //     Other, fancier options
00216 //    would be possible, at the cost of adding one more parameter to the PROFILE
00217 //    command. For example, if all Y variables are distributed according to some
00218 //    known Gaussian of standard deviation Sigma (which can be different for each measurement), 
00219 //    and the profile has been filled  with a weight equal to 1/Sigma**2, 
00220 //    then one cam use the following option: 
00221 // 
00222 //     'G'            Errors are 1./SQRT(Sum(1/sigma**2)) 
00223 //    For example, this would be useful when all Y's are experimental quantities
00224 //    measured with different precision Sigma_Y.
00225 //
00226 //
00227 
00228    SetErrorOption(option);
00229 
00230    fBinEntries.Set(fNcells);  //*-* create number of entries per bin array
00231 
00232    // TH1::Sumw2 create sum of square of weights array times y (fSumw2) . This is always created for a TProfile
00233    TH1::Sumw2();                   //*-* create sum of squares of weights array times y
00234    // TProfile::Sumw2 create sum of square of weight2 (fBinSumw2). This is needed only for profile filled with weights not 1
00235    if (fgDefaultSumw2) Sumw2();    // optionally create sum of squares of weights / bin 
00236 
00237    fYmin = ymin;
00238    fYmax = ymax;
00239    fScaling = kFALSE;
00240    fTsumwy = fTsumwy2 = 0;
00241 
00242 }
00243 
00244 //______________________________________________________________________________
00245 TProfile::TProfile(const TProfile &profile) : TH1D()
00246 {
00247    // Copy constructor.
00248 
00249    ((TProfile&)profile).Copy(*this);
00250 }
00251 
00252 
00253 //______________________________________________________________________________
00254 void TProfile::Add(TF1 *, Double_t, Option_t * )
00255 {
00256    // Performs the operation: this = this + c1*f1
00257 
00258    Error("Add","Function not implemented for TProfile");
00259    return;
00260 }
00261 
00262 
00263 //______________________________________________________________________________
00264 void TProfile::Add(const TH1 *h1, Double_t c1)
00265 {
00266    // Performs the operation: this = this + c1*h1
00267 
00268    if (!h1) {
00269       Error("Add","Attempt to add a non-existing profile");
00270       return;
00271    }
00272    if (!h1->InheritsFrom(TProfile::Class())) {
00273       Error("Add","Attempt to add a non-profile object");
00274       return;
00275    }
00276    
00277    TProfileHelper::Add(this, this, h1, 1, c1);
00278 }
00279 
00280 //______________________________________________________________________________
00281 void TProfile::Add(const TH1 *h1, const TH1 *h2, Double_t c1, Double_t c2)
00282 {
00283 //*-*-*-*-*Replace contents of this profile by the addition of h1 and h2*-*-*
00284 //*-*      =============================================================
00285 //
00286 //   this = c1*h1 + c2*h2
00287 //
00288 //   c1 and c2 are considered as weights applied to the two summed profiles. 
00289 //   The operation acts therefore like merging the two profiles with a weight c1 and c2
00290 //
00291 
00292    if (!h1 || !h2) {
00293       Error("Add","Attempt to add a non-existing profile");
00294       return;
00295    }
00296    if (!h1->InheritsFrom(TProfile::Class())) {
00297       Error("Add","Attempt to add a non-profile object");
00298       return;
00299    }
00300    if (!h2->InheritsFrom(TProfile::Class())) {
00301       Error("Add","Attempt to add a non-profile object");
00302       return;
00303    }
00304    TProfileHelper::Add(this, h1, h2, c1, c2);
00305 }
00306 
00307 
00308 //______________________________________________________________________________
00309 void TProfile::Approximate(Bool_t approx)
00310 {
00311 //     static function
00312 // set the fgApproximate flag. When the flag is true, the function GetBinError
00313 // will approximate the bin error with the average profile error on all bins
00314 // in the following situation only
00315 //  - the number of bins in the profile is less than 1002
00316 //  - the bin number of entries is small ( <5)
00317 //  - the estimated bin error is extremely small compared to the bin content
00318 //  (see TProfile::GetBinError)
00319 
00320    fgApproximate = approx;
00321 }
00322 
00323 //______________________________________________________________________________
00324 Int_t TProfile::BufferEmpty(Int_t action)
00325 {
00326 // Fill histogram with all entries in the buffer.
00327 // action = -1 histogram is reset and refilled from the buffer (called by THistPainter::Paint)
00328 // action =  0 histogram is filled from the buffer
00329 // action =  1 histogram is filled and buffer is deleted
00330 //             The buffer is automatically deleted when the number of entries
00331 //             in the buffer is greater than the number of entries in the histogram
00332 
00333    // do we need to compute the bin size?
00334    if (!fBuffer) return 0;
00335    Int_t nbentries = (Int_t)fBuffer[0];
00336    if (!nbentries) return 0;
00337    Double_t *buffer = fBuffer;
00338    if (nbentries < 0) {
00339       if (action == 0) return 0;
00340       nbentries  = -nbentries;
00341       fBuffer=0;
00342       Reset();
00343       fBuffer = buffer;
00344    }
00345    if (TestBit(kCanRebin) || fXaxis.GetXmax() <= fXaxis.GetXmin()) {
00346       //find min, max of entries in buffer
00347       Double_t xmin = fBuffer[2];
00348       Double_t xmax = xmin;
00349       for (Int_t i=1;i<nbentries;i++) {
00350          Double_t x = fBuffer[3*i+2];
00351          if (x < xmin) xmin = x;
00352          if (x > xmax) xmax = x;
00353       }
00354       if (fXaxis.GetXmax() <= fXaxis.GetXmin()) {
00355          THLimitsFinder::GetLimitsFinder()->FindGoodLimits(this,xmin,xmax);
00356       } else {
00357          fBuffer = 0;
00358          Int_t keep = fBufferSize; fBufferSize = 0;
00359          if (xmin <  fXaxis.GetXmin()) RebinAxis(xmin,&fXaxis);
00360          if (xmax >= fXaxis.GetXmax()) RebinAxis(xmax,&fXaxis);
00361          fBuffer = buffer;
00362          fBufferSize = keep;
00363       }
00364    }
00365 
00366    fBuffer = 0;
00367 
00368    for (Int_t i=0;i<nbentries;i++) {
00369       Fill(buffer[3*i+2],buffer[3*i+3],buffer[3*i+1]);
00370    }
00371    fBuffer = buffer;
00372 
00373    if (action > 0) { delete [] fBuffer; fBuffer = 0; fBufferSize = 0;}
00374    else {
00375       if (nbentries == (Int_t)fEntries) fBuffer[0] = -nbentries;
00376       else                              fBuffer[0] = 0;
00377    }
00378    return nbentries;
00379 }
00380 
00381 //______________________________________________________________________________
00382 Int_t TProfile::BufferFill(Double_t x, Double_t y, Double_t w)
00383 {
00384 // accumulate arguments in buffer. When buffer is full, empty the buffer
00385 // fBuffer[0] = number of entries in buffer
00386 // fBuffer[1] = w of first entry
00387 // fBuffer[2] = x of first entry
00388 // fBuffer[3] = y of first entry
00389 
00390    if (!fBuffer) return -2;
00391    Int_t nbentries = (Int_t)fBuffer[0];
00392    if (nbentries < 0) {
00393       nbentries  = -nbentries;
00394       fBuffer[0] =  nbentries;
00395       if (fEntries > 0) {
00396          Double_t *buffer = fBuffer; fBuffer=0;
00397          Reset();
00398          fBuffer = buffer;
00399       }
00400    }
00401    if (3*nbentries+3 >= fBufferSize) {
00402       BufferEmpty(1);
00403       return Fill(x,y,w);
00404    }
00405    fBuffer[3*nbentries+1] = w;
00406    fBuffer[3*nbentries+2] = x;
00407    fBuffer[3*nbentries+3] = y;
00408    fBuffer[0] += 1;
00409    return -2;
00410 }
00411 
00412 //______________________________________________________________________________
00413 void TProfile::Copy(TObject &obj) const
00414 {
00415 //*-*-*-*-*-*-*-*Copy a Profile histogram to a new profile histogram*-*-*-*-*
00416 //*-*            ===================================================
00417 
00418    TH1D::Copy(((TProfile&)obj));
00419    fBinEntries.Copy(((TProfile&)obj).fBinEntries);
00420    fBinSumw2.Copy(((TProfile&)obj).fBinSumw2);
00421    for (int bin=0;bin<fNcells;bin++) {
00422       ((TProfile&)obj).fArray[bin]        = fArray[bin];
00423       ((TProfile&)obj).fSumw2.fArray[bin] = fSumw2.fArray[bin];
00424    }
00425 
00426    ((TProfile&)obj).fYmin = fYmin;
00427    ((TProfile&)obj).fYmax = fYmax;
00428    ((TProfile&)obj).fScaling   = fScaling;
00429    ((TProfile&)obj).fErrorMode = fErrorMode;
00430    ((TProfile&)obj).fTsumwy    = fTsumwy;
00431    ((TProfile&)obj).fTsumwy2   = fTsumwy2;
00432 }
00433 
00434 
00435 //______________________________________________________________________________
00436 void TProfile::Divide(TF1 *, Double_t )
00437 {
00438    // Performs the operation: this = this/(c1*f1)
00439 
00440    Error("Divide","Function not implemented for TProfile");
00441    return;
00442 }
00443 
00444 //______________________________________________________________________________
00445 void TProfile::Divide(const TH1 *h1)
00446 {
00447 //*-*-*-*-*-*-*-*-*-*-*Divide this profile by h1*-*-*-*-*-*-*-*-*-*-*-*-*
00448 //*-*                  =========================
00449 //
00450 //   this = this/h1
00451 // This function accepts to divide a TProfile by a histogram
00452 //
00453 
00454    if (!h1) {
00455       Error("Divide","Attempt to divide a non-existing profile");
00456       return;
00457    }
00458    if (!h1->InheritsFrom(TH1::Class())) {
00459       Error("Divide","Attempt to divide by a non-profile or non-histogram object");
00460       return;
00461    }
00462    TProfile *p1 = (TProfile*)h1;
00463 
00464    Int_t nbinsx = GetNbinsX();
00465 //*-*- Check profile compatibility
00466    if (nbinsx != p1->GetNbinsX()) {
00467       Error("Divide","Attempt to divide profiles with different number of bins");
00468       return;
00469    }
00470 
00471 //*-*- Reset statistics
00472    fEntries = fTsumw   = fTsumw2 = fTsumwx = fTsumwx2 = fTsumwy = fTsumwy2 = 0;
00473 
00474 //*-*- Loop on bins (including underflows/overflows)
00475    Int_t bin;
00476    Double_t *cu1=0, *er1=0, *en1=0;
00477    Double_t e0,e1,c12;
00478    if (h1->InheritsFrom(TProfile::Class())) {
00479       cu1 = p1->GetW();
00480       er1 = p1->GetW2();
00481       en1 = p1->GetB();
00482    }
00483    Double_t c0,c1,w,z,x;
00484    for (bin=0;bin<=nbinsx+1;bin++) {
00485       c0  = fArray[bin];
00486       if (cu1) c1  = cu1[bin];
00487       else     c1  = h1->GetBinContent(bin);
00488       if (c1) w = c0/c1;
00489       else    w = 0;
00490       fArray[bin] = w;
00491       z = TMath::Abs(w);
00492       x = fXaxis.GetBinCenter(bin);
00493       fEntries++;
00494       fTsumw   += z;
00495       fTsumw2  += z*z;
00496       fTsumwx  += z*x;
00497       fTsumwx2 += z*x*x;
00498       fTsumwy  += z*c1;
00499       fTsumwx2 += z*c1*c1;
00500       e0 = fSumw2.fArray[bin];
00501       if (er1) e1 = er1[bin];
00502       else    {e1 = h1->GetBinError(bin); e1*=e1;}
00503       c12= c1*c1;
00504       if (!c1) fSumw2.fArray[bin] = 0;
00505       else     fSumw2.fArray[bin] = (e0*c1*c1 + e1*c0*c0)/(c12*c12);
00506       if (!en1) continue;
00507       if (!en1[bin]) fBinEntries.fArray[bin] = 0;
00508       else           fBinEntries.fArray[bin] /= en1[bin];
00509    }
00510    // mantaining the correct sum of weights square is not supported when dividing
00511    // bin error resulting from division of profile needs to be checked 
00512    if (fBinSumw2.fN) { 
00513       Warning("Divide","Cannot preserve during the division of profiles the sum of bin weight square");
00514       fBinSumw2 = TArrayD();
00515    }
00516    
00517 }
00518 
00519 
00520 //______________________________________________________________________________
00521 void TProfile::Divide(const TH1 *h1, const TH1 *h2, Double_t c1, Double_t c2, Option_t *option)
00522 {
00523 //*-*-*-*-*Replace contents of this profile by the division of h1 by h2*-*-*
00524 //*-*      ============================================================
00525 //
00526 //   this = c1*h1/(c2*h2)
00527 //
00528 
00529    TString opt = option;
00530    opt.ToLower();
00531    Bool_t binomial = kFALSE;
00532    if (opt.Contains("b")) binomial = kTRUE;
00533    if (!h1 || !h2) {
00534       Error("Divide","Attempt to divide a non-existing profile");
00535       return;
00536    }
00537    if (!h1->InheritsFrom(TProfile::Class())) {
00538       Error("Divide","Attempt to divide a non-profile object");
00539       return;
00540    }
00541    TProfile *p1 = (TProfile*)h1;
00542    if (!h2->InheritsFrom(TProfile::Class())) {
00543       Error("Divide","Attempt to divide by a non-profile object");
00544       return;
00545    }
00546    TProfile *p2 = (TProfile*)h2;
00547 
00548    Int_t nbinsx = GetNbinsX();
00549 //*-*- Check histogram compatibility
00550    if (nbinsx != p1->GetNbinsX() || nbinsx != p2->GetNbinsX()) {
00551       Error("Divide","Attempt to divide profiles with different number of bins");
00552       return;
00553    }
00554    if (!c2) {
00555       Error("Divide","Coefficient of dividing profile cannot be zero");
00556       return;
00557    }
00558 
00559    //THE ALGORITHM COMPUTING THE ERRORS IS WRONG. HELP REQUIRED
00560    printf("WARNING!!: The algorithm in TProfile::Divide computing the errors is not accurate\n");
00561    printf(" Instead of Divide(TProfile *h1, TProfile *h2, do:\n");
00562    printf("   TH1D *p1 = h1->ProjectionX();\n");
00563    printf("   TH1D *p2 = h2->ProjectionX();\n");
00564    printf("   p1->Divide(p2);\n");
00565 
00566 //*-*- Reset statistics
00567    fEntries = fTsumw   = fTsumw2 = fTsumwx = fTsumwx2 = 0;
00568 
00569 //*-*- Loop on bins (including underflows/overflows)
00570    Int_t bin;
00571    Double_t *cu1 = p1->GetW();
00572    Double_t *cu2 = p2->GetW();
00573    Double_t *er1 = p1->GetW2();
00574    Double_t *er2 = p2->GetW2();
00575    Double_t *en1 = p1->GetB();
00576    Double_t *en2 = p2->GetB();
00577    Double_t b1,b2,w,z,x,ac1,ac2;
00578    //d1 = c1*c1;
00579    //d2 = c2*c2;
00580    ac1 = TMath::Abs(c1);
00581    ac2 = TMath::Abs(c2);
00582    for (bin=0;bin<=nbinsx+1;bin++) {
00583       b1  = cu1[bin];
00584       b2  = cu2[bin];
00585       if (b2) w = c1*b1/(c2*b2);
00586       else    w = 0;
00587       fArray[bin] = w;
00588       z = TMath::Abs(w);
00589       x = fXaxis.GetBinCenter(bin);
00590       fEntries++;
00591       fTsumw   += z;
00592       fTsumw2  += z*z;
00593       fTsumwx  += z*x;
00594       fTsumwx2 += z*x*x;
00595       //fTsumwy  += z*x;
00596       //fTsumwy2 += z*x*x;
00597       Double_t e1 = er1[bin];
00598       Double_t e2 = er2[bin];
00599     //Double_t b22= b2*b2*d2;
00600       Double_t b22= b2*b2*TMath::Abs(c2);
00601       if (!b2) fSumw2.fArray[bin] = 0;
00602       else {
00603          if (binomial) {
00604             fSumw2.fArray[bin] = TMath::Abs(w*(1-w)/b2);
00605          } else {
00606           //fSumw2.fArray[bin] = d1*d2*(e1*b2*b2 + e2*b1*b1)/(b22*b22);
00607             fSumw2.fArray[bin] = ac1*ac2*(e1*b2*b2 + e2*b1*b1)/(b22*b22);
00608          }
00609       }
00610       if (en2[bin]) fBinEntries.fArray[bin] = en1[bin]/en2[bin];
00611       else          fBinEntries.fArray[bin] = 0;
00612    }
00613 
00614    // mantaining the correct sum of weights square is not supported when dividing
00615    // bin error resulting from division of profile needs to be checked 
00616    if (fBinSumw2.fN) { 
00617       Warning("Divide","Cannot preserve during the division of profiles the sum of bin weight square");
00618       fBinSumw2 = TArrayD();
00619    }
00620 
00621 }
00622 
00623 //______________________________________________________________________________
00624 TH1 *TProfile::DrawCopy(Option_t *option) const
00625 {
00626 //*-*-*-*-*-*-*-*Draw a copy of this profile histogram*-*-*-*-*-*-*-*-*-*-*-*
00627 //*-*            =====================================
00628    TString opt = option;
00629    opt.ToLower();
00630    if (gPad && !opt.Contains("same")) gPad->Clear();
00631    TProfile *newpf = (TProfile*)Clone();
00632    newpf->SetDirectory(0);
00633    newpf->SetBit(kCanDelete);
00634    newpf->AppendPad(option);
00635    return newpf;
00636 }
00637 
00638 //______________________________________________________________________________
00639 Int_t TProfile::Fill(Double_t x, Double_t y)
00640 {
00641 //*-*-*-*-*-*-*-*-*-*-*Fill a Profile histogram (no weights)*-*-*-*-*-*-*-*
00642 //*-*                  =====================================
00643 
00644    if (fBuffer) return BufferFill(x,y,1);
00645 
00646    Int_t bin;
00647    if (fYmin != fYmax) {
00648       if (y <fYmin || y> fYmax) return -1;
00649    }
00650 
00651    fEntries++;
00652    bin =fXaxis.FindBin(x);
00653    AddBinContent(bin, y);
00654    fSumw2.fArray[bin] += (Double_t)y*y;
00655    fBinEntries.fArray[bin] += 1;
00656    if (fBinSumw2.fN)  fBinSumw2.fArray[bin] += 1;
00657    if (bin == 0 || bin > fXaxis.GetNbins()) {
00658       if (!fgStatOverflows) return -1;
00659    }
00660    fTsumw++;
00661    fTsumw2++;
00662    fTsumwx  += x;
00663    fTsumwx2 += x*x;
00664    fTsumwy  += y;
00665    fTsumwy2 += y*y;
00666    return bin;
00667 }
00668 
00669 //______________________________________________________________________________
00670 Int_t TProfile::Fill(const char *namex, Double_t y)
00671 {
00672 // Fill a Profile histogram (no weights)
00673 //
00674    Int_t bin;
00675    if (fYmin != fYmax) {
00676       if (y <fYmin || y> fYmax) return -1;
00677    }
00678 
00679    fEntries++;
00680    bin =fXaxis.FindBin(namex);
00681    AddBinContent(bin, y);
00682    fSumw2.fArray[bin] += (Double_t)y*y;
00683    fBinEntries.fArray[bin] += 1;
00684    if (fBinSumw2.fN)  fBinSumw2.fArray[bin] += 1;
00685    if (bin == 0 || bin > fXaxis.GetNbins()) {
00686       if (!fgStatOverflows) return -1;
00687    }
00688    Double_t x = fXaxis.GetBinCenter(bin);
00689    fTsumw++;
00690    fTsumw2++;
00691    fTsumwx  += x;
00692    fTsumwx2 += x*x;
00693    fTsumwy  += y;
00694    fTsumwy2 += y*y;
00695    return bin;
00696 }
00697 
00698 //______________________________________________________________________________
00699 Int_t TProfile::Fill(Double_t x, Double_t y, Double_t w)
00700 {
00701 //*-*-*-*-*-*-*-*-*-*-*Fill a Profile histogram with weights*-*-*-*-*-*-*-*
00702 //*-*                  =====================================
00703 
00704    if (fBuffer) return BufferFill(x,y,w);
00705 
00706    Int_t bin;
00707    if (fYmin != fYmax) {
00708       if (y <fYmin || y> fYmax) return -1;
00709    }
00710 
00711    Double_t u= (w > 0 ? w : -w);
00712    fEntries++;
00713    bin =fXaxis.FindBin(x);
00714    AddBinContent(bin, u*y);
00715    fSumw2.fArray[bin] += u*y*y;
00716    fBinEntries.fArray[bin] += u;
00717    if (fBinSumw2.fN)  fBinSumw2.fArray[bin] += u*u;
00718    if (bin == 0 || bin > fXaxis.GetNbins()) {
00719       if (!fgStatOverflows) return -1;
00720    }
00721    fTsumw   += u;
00722    fTsumw2  += u*u;
00723    fTsumwx  += u*x;
00724    fTsumwx2 += u*x*x;
00725    fTsumwy  += u*y;
00726    fTsumwy2 += u*y*y;
00727    return bin;
00728 }
00729 
00730 //______________________________________________________________________________
00731 Int_t TProfile::Fill(const char *namex, Double_t y, Double_t w)
00732 {
00733 // Fill a Profile histogram with weights
00734 //
00735    Int_t bin;
00736 
00737    if (fYmin != fYmax) {
00738       if (y <fYmin || y> fYmax) return -1;
00739    }
00740 
00741    Double_t u= (w > 0 ? w : -w);
00742    fEntries++;
00743    bin =fXaxis.FindBin(namex);
00744    AddBinContent(bin, u*y);
00745    fSumw2.fArray[bin] += u*y*y;
00746    fBinEntries.fArray[bin] += u;
00747    if (fBinSumw2.fN)  fBinSumw2.fArray[bin] += u*u;
00748    if (bin == 0 || bin > fXaxis.GetNbins()) {
00749       if (!fgStatOverflows) return -1;
00750    }
00751    Double_t x = fXaxis.GetBinCenter(bin);
00752    fTsumw   += u;
00753    fTsumw2  += u*u;
00754    fTsumwx  += u*x;
00755    fTsumwx2 += u*x*x;
00756    fTsumwy  += u*y;
00757    fTsumwy2 += u*y*y;
00758    return bin;
00759 }
00760 
00761 
00762 //______________________________________________________________________________
00763 void TProfile::FillN(Int_t ntimes, const Double_t *x, const Double_t *y, const Double_t *w, Int_t stride)
00764 {
00765 //*-*-*-*-*-*-*-*-*-*-*Fill a Profile histogram with weights*-*-*-*-*-*-*-*
00766 //*-*                  =====================================
00767    Int_t bin,i;
00768    ntimes *= stride;
00769    for (i=0;i<ntimes;i+=stride) {
00770       if (fYmin != fYmax) {
00771          if (y[i] <fYmin || y[i]> fYmax) continue;
00772       }
00773 
00774       Double_t u= (w[i] > 0 ? w[i] : -w[i]);
00775       fEntries++;
00776       bin =fXaxis.FindBin(x[i]);
00777       AddBinContent(bin, u*y[i]);
00778       fSumw2.fArray[bin] += u*y[i]*y[i];
00779       fBinEntries.fArray[bin] += u;
00780       if (fBinSumw2.fN)  fBinSumw2.fArray[bin] += u*u;
00781       if (bin == 0 || bin > fXaxis.GetNbins()) {
00782          if (!fgStatOverflows) continue;
00783       }
00784       fTsumw   += u;
00785       fTsumw2  += u*u;
00786       fTsumwx  += u*x[i];
00787       fTsumwx2 += u*x[i]*x[i];
00788       fTsumwy  += u*y[i];
00789       fTsumwy2 += u*y[i]*y[i];
00790    }
00791 }
00792 
00793 //______________________________________________________________________________
00794 Double_t TProfile::GetBinContent(Int_t bin) const
00795 {
00796 //*-*-*-*-*-*-*Return bin content of a Profile histogram*-*-*-*-*-*-*-*-*-*
00797 //*-*          =========================================
00798 
00799    if (fBuffer) ((TProfile*)this)->BufferEmpty();
00800 
00801    if (bin < 0 || bin >= fNcells) return 0;
00802    if (fBinEntries.fArray[bin] == 0) return 0;
00803    if (!fArray) return 0;
00804    return fArray[bin]/fBinEntries.fArray[bin];
00805 }
00806 
00807 //______________________________________________________________________________
00808 Double_t TProfile::GetBinEntries(Int_t bin) const
00809 {
00810 //*-*-*-*-*-*-*Return bin entries of a Profile histogram*-*-*-*-*-*-*-*-*-*
00811 //*-*          =========================================
00812 
00813    if (fBuffer) ((TProfile*)this)->BufferEmpty();
00814 
00815    if (bin < 0 || bin >= fNcells) return 0;
00816    return fBinEntries.fArray[bin];
00817 }
00818 
00819 //______________________________________________________________________________
00820 Double_t TProfile::GetBinEffectiveEntries(Int_t bin) const
00821 {
00822 //            Return bin effective entries for a weighted filled Profile histogram. 
00823 //            In case of an unweighted profile, it is equivalent to the number of entries per bin   
00824 //            The effective entries is defined as the square of the sum of the weights divided by the 
00825 //            sum of the weights square. 
00826 //            TProfile::Sumw2() must be called before filling the profile with weights. 
00827 //            Only by calling this method the  sum of the square of the weights per bin is stored. 
00828 //  
00829 //*-*          =========================================
00830 
00831    return TProfileHelper::GetBinEffectiveEntries((TProfile*)this, bin);
00832 }
00833 
00834 //______________________________________________________________________________
00835 Double_t TProfile::GetBinError(Int_t bin) const
00836 {
00837 // *-*-*-*-*-*-*Return bin error of a Profile histogram*-*-*-*-*-*-*-*-*-*
00838 // *-*          =======================================
00839 //
00840 // Computing errors: A moving field
00841 // =================================
00842 // The computation of errors for a TProfile has evolved with the versions
00843 // of ROOT. The difficulty is in computing errors for bins with low statistics.
00844 // - prior to version 3.00, we had no special treatment of low statistic bins.
00845 //   As a result, these bins had huge errors. The reason is that the
00846 //   expression eprim2 is very close to 0 (rounding problems) or 0.
00847 // - in version 3.00 (18 Dec 2000), the algorithm is protected for values of
00848 //   eprim2 very small and the bin errors set to the average bin errors, following
00849 //   recommendations from a group of users.
00850 // - in version 3.01 (19 Apr 2001), it is realized that the algorithm above
00851 //   should be applied only to low statistic bins.
00852 // - in version 3.02 (26 Sep 2001), the same group of users recommend instead
00853 //   to take two times the average error on all bins for these low
00854 //   statistics bins giving a very small value for eprim2.
00855 // - in version 3.04 (Nov 2002), the algorithm is modified/protected for the case
00856 //   when a TProfile is projected (ProjectionX). The previous algorithm
00857 //   generated a N^2 problem when projecting a TProfile with a large number of
00858 //   bins (eg 100000).
00859 // - in version 3.05/06, a new static function TProfile::Approximate
00860 //   is introduced to enable or disable (default) the approximation.
00861 //
00862 // Ideas for improvements of this algorithm are welcome. No suggestions
00863 // received since our call for advice to roottalk in Jul 2002.
00864 // see for instance: http://root.cern.ch/root/roottalk/roottalk02/2916.html
00865 
00866    return TProfileHelper::GetBinError((TProfile*)this, bin);
00867 }
00868 
00869 //______________________________________________________________________________
00870 Option_t *TProfile::GetErrorOption() const
00871 {
00872 //*-*-*-*-*-*-*-*-*-*Return option to compute profile errors*-*-*-*-*-*-*-*-*
00873 //*-*                =======================================
00874 
00875    if (fErrorMode == kERRORSPREAD)  return "s";
00876    if (fErrorMode == kERRORSPREADI) return "i";
00877    if (fErrorMode == kERRORSPREADG) return "g";
00878    return "";
00879 }
00880 
00881 //______________________________________________________________________________
00882 char* TProfile::GetObjectInfo(Int_t px, Int_t py) const
00883 {
00884    //   Redefines TObject::GetObjectInfo.
00885    //   Displays the profile info (bin number, contents, eroor, entries per bin
00886    //   corresponding to cursor position px,py
00887    //
00888    if (!gPad) return (char*)"";
00889    static char info[200];
00890    Double_t x  = gPad->PadtoX(gPad->AbsPixeltoX(px));
00891    Double_t y  = gPad->PadtoY(gPad->AbsPixeltoY(py));
00892    Int_t binx   = GetXaxis()->FindFixBin(x);
00893    snprintf(info,200,"(x=%g, y=%g, binx=%d, binc=%g, bine=%g, binn=%d)", x, y, binx, GetBinContent(binx), GetBinError(binx), (Int_t)GetBinEntries(binx));
00894    return info;
00895 }
00896 
00897 //______________________________________________________________________________
00898 void TProfile::GetStats(Double_t *stats) const
00899 {
00900    // fill the array stats from the contents of this profile
00901    // The array stats must be correctly dimensionned in the calling program.
00902    // stats[0] = sumw
00903    // stats[1] = sumw2
00904    // stats[2] = sumwx
00905    // stats[3] = sumwx2
00906    // stats[4] = sumwy
00907    // stats[5] = sumwy2
00908    //
00909    // If no axis-subrange is specified (via TAxis::SetRange), the array stats
00910    // is simply a copy of the statistics quantities computed at filling time.
00911    // If a sub-range is specified, the function recomputes these quantities
00912    // from the bin contents in the current axis range.
00913 
00914    if (fBuffer) ((TProfile*)this)->BufferEmpty();
00915 
00916    // Loop on bins
00917    Int_t bin, binx;
00918    if (fTsumw == 0 || fXaxis.TestBit(TAxis::kAxisRange)) {
00919       for (bin=0;bin<6;bin++) stats[bin] = 0;
00920       if (!fBinEntries.fArray) return;
00921       Int_t firstBinX = fXaxis.GetFirst();
00922       Int_t lastBinX  = fXaxis.GetLast();
00923       // include underflow/overflow if TH1::StatOverflows(kTRUE) in case no range is set on the axis
00924       if (fgStatOverflows && !fXaxis.TestBit(TAxis::kAxisRange)) {
00925          if (firstBinX == 1) firstBinX = 0;
00926          if (lastBinX ==  fXaxis.GetNbins() ) lastBinX += 1;
00927       }
00928       for (binx = firstBinX; binx <= lastBinX; binx++) {
00929          Double_t w   = fBinEntries.fArray[binx];
00930          Double_t w2  = (fBinSumw2.fN ? fBinSumw2.fArray[binx] : w);  
00931          Double_t x   = fXaxis.GetBinCenter(binx);
00932          stats[0] += w;
00933          stats[1] += w2;
00934          stats[2] += w*x;
00935          stats[3] += w*x*x;
00936          stats[4] += fArray[binx];
00937          stats[5] += fSumw2.fArray[binx];
00938       }
00939    } else {
00940       if (fTsumwy == 0 && fTsumwy2 == 0) {
00941          //this case may happen when processing TProfiles with version <=3
00942          TProfile *p = (TProfile*)this; // chheting with const
00943          for (binx=fXaxis.GetFirst();binx<=fXaxis.GetLast();binx++) {
00944             p->fTsumwy  += fArray[binx];
00945             p->fTsumwy2 += fSumw2.fArray[binx];
00946          }
00947       }
00948       stats[0] = fTsumw;
00949       stats[1] = fTsumw2;
00950       stats[2] = fTsumwx;
00951       stats[3] = fTsumwx2;
00952       stats[4] = fTsumwy;
00953       stats[5] = fTsumwy2;
00954    }
00955 }
00956 
00957 //___________________________________________________________________________
00958 void TProfile::LabelsDeflate(Option_t *option)
00959 {
00960 // Reduce the number of bins for this axis to the number of bins having a label.
00961 
00962    TProfileHelper::LabelsDeflate(this, option);
00963 }
00964 
00965 //___________________________________________________________________________
00966 void TProfile::LabelsInflate(Option_t *options)
00967 {
00968 // Double the number of bins for axis.
00969 // Refill histogram
00970 // This function is called by TAxis::FindBin(const char *label)
00971 
00972   TProfileHelper::LabelsInflate(this, options); 
00973 }
00974 
00975 //___________________________________________________________________________
00976 void TProfile::LabelsOption(Option_t *option, Option_t * /*ax */)
00977 {
00978 //  Set option(s) to draw axis with labels
00979 //  option = "a" sort by alphabetic order
00980 //         = ">" sort by decreasing values
00981 //         = "<" sort by increasing values
00982 //         = "h" draw labels horizonthal
00983 //         = "v" draw labels vertical
00984 //         = "u" draw labels up (end of label right adjusted)
00985 //         = "d" draw labels down (start of label left adjusted)
00986 
00987    THashList *labels = fXaxis.GetLabels();
00988    if (!labels) {
00989       Warning("LabelsOption","Cannot sort. No labels");
00990       return;
00991    }
00992    TString opt = option;
00993    opt.ToLower();
00994    if (opt.Contains("h")) {
00995       fXaxis.SetBit(TAxis::kLabelsHori);
00996       fXaxis.ResetBit(TAxis::kLabelsVert);
00997       fXaxis.ResetBit(TAxis::kLabelsDown);
00998       fXaxis.ResetBit(TAxis::kLabelsUp);
00999    }
01000    if (opt.Contains("v")) {
01001       fXaxis.SetBit(TAxis::kLabelsVert);
01002       fXaxis.ResetBit(TAxis::kLabelsHori);
01003       fXaxis.ResetBit(TAxis::kLabelsDown);
01004       fXaxis.ResetBit(TAxis::kLabelsUp);
01005    }
01006    if (opt.Contains("u")) {
01007       fXaxis.SetBit(TAxis::kLabelsUp);
01008       fXaxis.ResetBit(TAxis::kLabelsVert);
01009       fXaxis.ResetBit(TAxis::kLabelsDown);
01010       fXaxis.ResetBit(TAxis::kLabelsHori);
01011    }
01012    if (opt.Contains("d")) {
01013       fXaxis.SetBit(TAxis::kLabelsDown);
01014       fXaxis.ResetBit(TAxis::kLabelsVert);
01015       fXaxis.ResetBit(TAxis::kLabelsHori);
01016       fXaxis.ResetBit(TAxis::kLabelsUp);
01017    }
01018    Int_t sort = -1;
01019    if (opt.Contains("a")) sort = 0;
01020    if (opt.Contains(">")) sort = 1;
01021    if (opt.Contains("<")) sort = 2;
01022    if (sort < 0) return;
01023 
01024    Int_t n = TMath::Min(fXaxis.GetNbins(), labels->GetSize());
01025    Int_t *a = new Int_t[n+2];
01026    Int_t i,j;
01027    Double_t *cont   = new Double_t[n+2];
01028    Double_t *sumw   = new Double_t[n+2];
01029    Double_t *errors = new Double_t[n+2];
01030    Double_t *ent    = new Double_t[n+2];
01031    THashList *labold = new THashList(labels->GetSize(),1);
01032    TIter nextold(labels);
01033    TObject *obj;
01034    while ((obj=nextold())) {
01035       labold->Add(obj);
01036    }
01037    labels->Clear();
01038    if (sort > 0) {
01039       //---sort by values of bins
01040       for (i=1;i<=n;i++) {
01041          sumw[i-1]   = fArray[i];
01042          errors[i-1] = fSumw2.fArray[i];
01043          ent[i-1]    = fBinEntries.fArray[i];
01044          if (fBinEntries.fArray[i] == 0) cont[i-1] = 0;
01045          else cont[i-1] = fArray[i]/fBinEntries.fArray[i];
01046       }
01047       if (sort ==1) TMath::Sort(n,cont,a,kTRUE);  //sort by decreasing values
01048       else          TMath::Sort(n,cont,a,kFALSE); //sort by increasing values
01049       for (i=1;i<=n;i++) {
01050          fArray[i] = sumw[a[i-1]];
01051          fSumw2.fArray[i] = errors[a[i-1]];
01052          fBinEntries.fArray[i] = ent[a[i-1]];
01053       }
01054       for (i=1;i<=n;i++) {
01055          obj = labold->At(a[i-1]);
01056          labels->Add(obj);
01057          obj->SetUniqueID(i);
01058       }
01059    } else {
01060       //---alphabetic sort
01061       const UInt_t kUsed = 1<<18;
01062       TObject *objk=0;
01063       a[0] = 0;
01064       a[n+1] = n+1;
01065       for (i=1;i<=n;i++) {
01066          const char *label = "zzzzzzzzzzzz";
01067          for (j=1;j<=n;j++) {
01068             obj = labold->At(j-1);
01069             if (!obj) continue;
01070             if (obj->TestBit(kUsed)) continue;
01071             //use strcasecmp for case non-sensitive sort (may be an option)
01072             if (strcmp(label,obj->GetName()) < 0) continue;
01073             objk = obj;
01074             a[i] = j;
01075             label = obj->GetName();
01076          }
01077          if (objk) {
01078             objk->SetUniqueID(i);
01079             labels->Add(objk);
01080             objk->SetBit(kUsed);
01081          }
01082       }
01083       for (i=1;i<=n;i++) {
01084          obj = labels->At(i-1);
01085          if (!obj) continue;
01086          obj->ResetBit(kUsed);
01087       }
01088 
01089       for (i=1;i<=n;i++) {
01090          sumw[i]   = fArray[a[i]];
01091          errors[i] = fSumw2.fArray[a[i]];
01092          ent[i]    = fBinEntries.fArray[a[i]];
01093       }
01094       for (i=1;i<=n;i++) {
01095          fArray[i] = sumw[i];
01096          fSumw2.fArray[i] = errors[i];
01097          fBinEntries.fArray[i] = ent[i];
01098       }
01099    }
01100    delete labold;
01101    if (a)      delete [] a;
01102    if (sumw)   delete [] sumw;
01103    if (cont)   delete [] cont;
01104    if (errors) delete [] errors;
01105    if (ent)    delete [] ent;
01106 }
01107 
01108 //______________________________________________________________________________
01109 Long64_t TProfile::Merge(TCollection *li)
01110 {
01111    //Merge all histograms in the collection in this histogram.
01112    //This function computes the min/max for the x axis,
01113    //compute a new number of bins, if necessary,
01114    //add bin contents, errors and statistics.
01115    //If overflows are present and limits are different the function will fail.
01116    //The function returns the total number of entries in the result histogram
01117    //if the merge is successfull, -1 otherwise.
01118    //
01119    //IMPORTANT remark. The axis x may have different number
01120    //of bins and different limits, BUT the largest bin width must be
01121    //a multiple of the smallest bin width and the upper limit must also
01122    //be a multiple of the bin width.
01123 
01124    return TProfileHelper::Merge(this, li);
01125 }
01126 
01127 
01128 //______________________________________________________________________________
01129 void TProfile::Multiply(TF1 *f1, Double_t c1)
01130 {
01131    // Performs the operation: this = this*c1*f1
01132 
01133    if (!f1) {
01134       Error("Multiply","Attempt to multiply by a null function");
01135       return;
01136    }
01137 
01138    Int_t nbinsx = GetNbinsX();
01139 
01140 //*-*- Add statistics
01141    Double_t xx[1], cf1, ac1 = TMath::Abs(c1);
01142    Double_t s1[10];
01143    Int_t i;
01144    for (i=0;i<10;i++) {s1[i] = 0;}
01145    PutStats(s1);
01146 
01147    SetMinimum();
01148    SetMaximum();
01149 
01150 //*-*- Loop on bins (including underflows/overflows)
01151    Int_t bin;
01152    for (bin=0;bin<=nbinsx+1;bin++) {
01153       xx[0] = fXaxis.GetBinCenter(bin);
01154       if (!f1->IsInside(xx)) continue;
01155       TF1::RejectPoint(kFALSE);
01156       cf1 = f1->EvalPar(xx);
01157       if (TF1::RejectedPoint()) continue;
01158       fArray[bin]             *= c1*cf1;
01159       //see http://savannah.cern.ch/bugs/?func=detailitem&item_id=14851
01160       //fSumw2.fArray[bin]      *= c1*c1*cf1*cf1;
01161       fSumw2.fArray[bin]      *= ac1*cf1*cf1;
01162       //fBinEntries.fArray[bin] *= ac1*TMath::Abs(cf1);
01163    }
01164 }
01165 
01166 //______________________________________________________________________________
01167 void TProfile::Multiply(const TH1 *)
01168 {
01169 //*-*-*-*-*-*-*-*-*-*-*Multiply this profile by h1*-*-*-*-*-*-*-*-*-*-*-*-*
01170 //*-*                  =============================
01171 //
01172 //   this = this*h1
01173 //
01174    Error("Multiply","Multiplication of profile histograms not implemented");
01175 }
01176 
01177 
01178 //______________________________________________________________________________
01179 void TProfile::Multiply(const TH1 *, const TH1 *, Double_t, Double_t, Option_t *)
01180 {
01181 //*-*-*-*-*Replace contents of this profile by multiplication of h1 by h2*-*
01182 //*-*      ================================================================
01183 //
01184 //   this = (c1*h1)*(c2*h2)
01185 //
01186 
01187    Error("Multiply","Multiplication of profile histograms not implemented");
01188 }
01189 
01190 //______________________________________________________________________________
01191 TH1D *TProfile::ProjectionX(const char *name, Option_t *option) const
01192 {
01193 //*-*-*-*-*Project this profile into a 1-D histogram along X*-*-*-*-*-*-*
01194 //*-*      =================================================
01195 //
01196 //   The projection is always of the type TH1D.
01197 //
01198 //   if option "E" is specified the errors of the projected histogram are computed and set 
01199 //      to be equal to the errors of the profile.
01200 //      Option "E" is defined as the default one in the header file. 
01201 //   if option "" is specified the histogram errors are simply the sqrt of its content
01202 //   if option "B" is specified, the content of bin of the returned histogram
01203 //      will be equal to the GetBinEntries(bin) of the profile,
01204 //      otherwise (default) it will be equal to GetBinContent(bin)
01205 //   if option "C=E" the bin contents of the projection are set to the
01206 //       bin errors of the profile
01207 //   if option "W" is specified the bin content of the projected histogram  is set to the 
01208 //       product of the bin content of the profile and the entries. 
01209 //       With this option the returned histogram will be equivalent to the one obtained by 
01210 //       filling directly a TH1D using the 2-nd value as a weight. 
01211 //       This makes sense only for profile filled with weights =1. If not, the error of the 
01212 //        projected histogram obtained with this option will not be correct.
01213 
01214 
01215    TString opt = option;
01216    opt.ToLower();
01217    Int_t nx = fXaxis.GetNbins();
01218 
01219 // Create the projection histogram
01220    TString pname = name; 
01221    if (pname == "_px") { 
01222       pname = GetName(); 
01223       pname.Append("_px");
01224    }
01225    TH1D *h1;
01226    const TArrayD *bins = fXaxis.GetXbins();
01227    if (bins->fN == 0) {
01228       h1 = new TH1D(pname,GetTitle(),nx,fXaxis.GetXmin(),fXaxis.GetXmax());
01229    } else {
01230       h1 = new TH1D(pname,GetTitle(),nx,bins->fArray);
01231    }
01232    Bool_t computeErrors = kFALSE;
01233    Bool_t cequalErrors  = kFALSE;
01234    Bool_t binEntries    = kFALSE;
01235    Bool_t binWeight     = kFALSE;
01236    if (opt.Contains("b")) binEntries = kTRUE;
01237    if (opt.Contains("e")) computeErrors = kTRUE;
01238    if (opt.Contains("w")) binWeight = kTRUE;
01239    if (opt.Contains("c=e")) {cequalErrors = kTRUE; computeErrors=kFALSE;}
01240    if (computeErrors || binWeight ) h1->Sumw2();
01241 
01242    // Fill the projected histogram
01243    Double_t cont;
01244    for (Int_t bin =0;bin<=nx+1;bin++) {
01245 
01246       if (binEntries)         cont = GetBinEntries(bin);
01247       else if (cequalErrors)  cont = GetBinError(bin);
01248       else if (binWeight)     cont = fArray[bin];  // bin content * bin entries
01249       else                    cont = GetBinContent(bin);    // default case
01250       
01251       h1->SetBinContent(bin ,cont);
01252 
01253       // if option E projected histogram errors are same as profile
01254       if (computeErrors ) h1->SetBinError(bin , GetBinError(bin) );
01255       // in case of option W bin error is deduced from bin sum of z**2 values of profile
01256       // this is correct only if the profile is filled with weights =1
01257       if (binWeight) h1->SetBinError(bin , TMath::Sqrt(fSumw2.fArray[bin] ) );
01258       // in case of bin entries and h1 has sumw2 set, we need to set also the bin error
01259       if (binEntries && h1->GetSumw2() ) {
01260          Double_t err2;
01261          if (fBinSumw2.fN) 
01262             err2 = fBinSumw2.fArray[bin]; 
01263          else 
01264             err2 = cont; // this is fBinEntries.fArray[bin]
01265          h1->SetBinError(bin, TMath::Sqrt(err2 ) ); 
01266       }
01267 
01268    }
01269 
01270    // Copy the axis attributes and the axis labels if needed.
01271    h1->GetXaxis()->ImportAttributes(this->GetXaxis());
01272    h1->GetYaxis()->ImportAttributes(this->GetYaxis());
01273    THashList* labels=this->GetXaxis()->GetLabels();
01274    if (labels) {
01275       TIter iL(labels);
01276       TObjString* lb;
01277       Int_t i = 1;
01278       while ((lb=(TObjString*)iL())) {
01279          h1->GetXaxis()->SetBinLabel(i,lb->String().Data());
01280          i++;
01281       }
01282    }
01283    
01284    h1->SetEntries(fEntries);
01285    return h1;
01286 }
01287 
01288 //______________________________________________________________________________
01289 void TProfile::PutStats(Double_t *stats)
01290 {
01291    // Replace current statistics with the values in array stats
01292 
01293    fTsumw   = stats[0];
01294    fTsumw2  = stats[1];
01295    fTsumwx  = stats[2];
01296    fTsumwx2 = stats[3];
01297    fTsumwy  = stats[4];
01298    fTsumwy2 = stats[5];
01299 }
01300 
01301 //______________________________________________________________________________
01302 TH1 *TProfile::Rebin(Int_t ngroup, const char*newname, const Double_t *xbins)
01303 {
01304 //*-*-*-*-*Rebin this profile grouping ngroup bins together*-*-*-*-*-*-*-*-*
01305 //*-*      ================================================
01306 //  -case 1  xbins=0
01307 //   if newname is not blank a new temporary profile hnew is created.
01308 //   else the current profile is modified (default)
01309 //   The parameter ngroup indicates how many bins of this have to me merged
01310 //   into one bin of hnew
01311 //   If the original profile has errors stored (via Sumw2), the resulting
01312 //   profile has new errors correctly calculated.
01313 //
01314 //   examples: if hp is an existing TProfile histogram with 100 bins
01315 //     hp->Rebin();  //merges two bins in one in hp: previous contents of hp are lost
01316 //     hp->Rebin(5); //merges five bins in one in hp
01317 //     TProfile *hnew = hp->Rebin(5,"hnew"); // creates a new profile hnew
01318 //                                       //merging 5 bins of hp in one bin
01319 //
01320 //   NOTE:  If ngroup is not an exact divider of the number of bins,
01321 //          the top limit of the rebinned profile is changed
01322 //          to the upper edge of the bin=newbins*ngroup and the corresponding
01323 //          bins are added to the overflow bin.
01324 //          Statistics will be recomputed from the new bin contents.
01325 //
01326 //  -case 2  xbins!=0
01327 //   a new profile is created (you should specify newname).
01328 //   The parameter ngroup is the number of variable size bins in the created profile
01329 //   The array xbins must contain ngroup+1 elements that represent the low-edge
01330 //   of the bins.
01331 //   The data of the old bins are added to the new bin which contains the bin center
01332 //   of the old bins. It is possible that information from the old binning are attached
01333 //   to the under-/overflow bins of the new binning.
01334 //
01335 //   examples: if hp is an existing TProfile with 100 bins
01336 //     Double_t xbins[25] = {...} array of low-edges (xbins[25] is the upper edge of last bin
01337 //     hp->Rebin(24,"hpnew",xbins);  //creates a new variable bin size profile hpnew
01338 
01339    Int_t nbins    = fXaxis.GetNbins();
01340    Double_t xmin  = fXaxis.GetXmin();
01341    Double_t xmax  = fXaxis.GetXmax();
01342    if ((ngroup <= 0) || (ngroup > nbins)) {
01343       Error("Rebin", "Illegal value of ngroup=%d",ngroup);
01344       return 0;
01345    }
01346    if (!newname && xbins) {
01347       Error("Rebin","if xbins is specified, newname must be given");
01348       return 0;
01349    }
01350 
01351    Int_t newbins = nbins/ngroup;
01352    if (!xbins) { 
01353       Int_t nbg = nbins/ngroup;
01354       if (nbg*ngroup != nbins) {
01355          Warning("Rebin", "ngroup=%d must be an exact divider of nbins=%d",ngroup,nbins);
01356       }
01357    }
01358    else {   
01359    // in the case of xbins given (rebinning in variable bins) ngroup is the new number of bins.
01360    //  and number of grouped bins is not constant. 
01361    // when looping for setting the contents for the new histogram we 
01362    // need to loop on all bins of original histogram. Set then ngroup=nbins
01363       newbins = ngroup;
01364       ngroup = nbins;
01365    }
01366 
01367    // Save old bin contents into a new array
01368    Double_t *oldBins   = new Double_t[nbins+2];
01369    Double_t *oldCount  = new Double_t[nbins+2];
01370    Double_t *oldErrors = new Double_t[nbins+2];
01371    Double_t *oldBinw2  = (fBinSumw2.fN ? new Double_t[nbins+2] : 0  ); 
01372    Int_t bin, i;
01373    Double_t *cu1 = GetW();
01374    Double_t *er1 = GetW2();
01375    Double_t *en1 = GetB();
01376    Double_t *ew1 = GetB2();
01377    
01378    for (bin=0;bin<=nbins+1;bin++) {
01379       oldBins[bin]   = cu1[bin];
01380       oldCount[bin]  = en1[bin];
01381       oldErrors[bin] = er1[bin];
01382       if (ew1 && fBinSumw2.fN) oldBinw2[bin]  = ew1[bin];
01383    }
01384 
01385    // create a clone of the old histogram if newname is specified
01386    TProfile *hnew = this;
01387    if ((newname && strlen(newname) > 0) || xbins) {
01388       hnew = (TProfile*)Clone(newname);
01389    }
01390 
01391    // in case of ngroup not an excat divider of nbins, 
01392    // top limit is changed (see NOTE in method comment) 
01393    if(!xbins && (newbins*ngroup != nbins)) {
01394       xmax = fXaxis.GetBinUpEdge(newbins*ngroup);
01395       hnew->fTsumw = 0; //stats must be reset because top bins will be moved to overflow bin
01396    }
01397 
01398    // set correctly the axis and resizes the bin arrays
01399    if(!xbins && (fXaxis.GetXbins()->GetSize() > 0)){ 
01400       // for rebinning of variable bins in a constant group 
01401       Double_t *bins = new Double_t[newbins+1];
01402       for(i = 0; i <= newbins; ++i) bins[i] = fXaxis.GetBinLowEdge(1+i*ngroup);
01403       hnew->SetBins(newbins,bins); //this also changes the bin array's
01404       delete [] bins;
01405    } else if (xbins) { 
01406       // when rebinning in variable bins
01407       hnew->SetBins(newbins,xbins);
01408    } else {
01409       hnew->SetBins(newbins,xmin,xmax);
01410    }
01411 
01412    // merge bin contents ignoring now underflow/overflows
01413    if (fBinSumw2.fN) hnew->Sumw2();
01414 
01415    // Start merging only once the new lowest edge is reached
01416    Int_t startbin = 1;
01417    const Double_t newxmin = hnew->GetXaxis()->GetBinLowEdge(1);
01418    while( fXaxis.GetBinCenter(startbin) < newxmin && startbin <= nbins ) {
01419       startbin++;
01420    }
01421    
01422    Double_t *cu2 = hnew->GetW();
01423    Double_t *er2 = hnew->GetW2();
01424    Double_t *en2 = hnew->GetB();
01425    Double_t *ew2 = hnew->GetB2();
01426    Int_t oldbin = startbin;
01427    Double_t binContent, binCount, binError, binSumw2;
01428    for (bin = 1;bin<=newbins;bin++) {
01429       binContent = 0;
01430       binCount   = 0;
01431       binError   = 0;
01432       binSumw2   = 0;
01433 
01434       //for xbins != 0: ngroup == nbins
01435       Int_t imax = ngroup;
01436       Double_t xbinmax = hnew->GetXaxis()->GetBinUpEdge(bin);
01437       for (i=0;i<ngroup;i++) {
01438          if((hnew == this && (oldbin+i > nbins)) ||
01439             (hnew != this && (fXaxis.GetBinCenter(oldbin+i) > xbinmax)))
01440          {
01441             imax = i;
01442             break;
01443          }
01444 
01445          binContent += oldBins[oldbin+i];
01446          binCount   += oldCount[oldbin+i];
01447          binError   += oldErrors[oldbin+i];
01448          if (fBinSumw2.fN) binSumw2 += oldBinw2[oldbin+i];
01449       }
01450    
01451       cu2[bin] = binContent;
01452       er2[bin] = binError;
01453       en2[bin] = binCount;
01454       if (fBinSumw2.fN) ew2[bin] = binSumw2;
01455       oldbin += imax;
01456    }
01457    // set bin statistics for underflow bin
01458    binContent = 0;
01459    binCount   = 0;
01460    binError   = 0;
01461    binSumw2   = 0;
01462    for(i=0;i<startbin;i++)
01463    {
01464       binContent += oldBins[i];
01465       binCount   += oldCount[i];
01466       binError   += oldErrors[i];
01467       if (fBinSumw2.fN) binSumw2 += oldBinw2[i];
01468    }
01469    hnew->fArray[0] = binContent;
01470    hnew->fBinEntries[0] = binCount;
01471    hnew->fSumw2[0] = binError;
01472    if ( fBinSumw2.fN ) hnew->fBinSumw2[0] = binSumw2;
01473 
01474    // set bin statistics for overflow bin
01475    binContent = 0;
01476    binCount   = 0;
01477    binError   = 0;
01478    binSumw2   = 0;
01479    for(i=oldbin;i<=nbins+1;i++)
01480    {
01481       binContent += oldBins[i];
01482       binCount   += oldCount[i];
01483       binError   += oldErrors[i];
01484       if (fBinSumw2.fN) binSumw2 += oldBinw2[i];
01485    }
01486    hnew->fArray[newbins+1] = binContent;
01487    hnew->fBinEntries[newbins+1] = binCount;   
01488    hnew->fSumw2[newbins+1] = binError;
01489    if ( fBinSumw2.fN ) hnew->fBinSumw2[newbins+1] = binSumw2;
01490 
01491 
01492    delete [] oldBins;
01493    delete [] oldCount;
01494    delete [] oldErrors;
01495    if (oldBinw2) delete [] oldBinw2; 
01496    return hnew;
01497 }
01498 
01499 //______________________________________________________________________________
01500 void TProfile::RebinAxis(Double_t x, TAxis *axis)
01501 {
01502 // Profile histogram is resized along x axis such that x is in the axis range.
01503 // The new axis limits are recomputed by doubling iteratively
01504 // the current axis range until the specified value x is within the limits.
01505 // The algorithm makes a copy of the histogram, then loops on all bins
01506 // of the old histogram to fill the rebinned histogram.
01507 // Takes into account errors (Sumw2) if any.
01508 // The bit kCanRebin must be set before invoking this function.
01509 //  Ex:  h->SetBit(TH1::kCanRebin);
01510 
01511    TProfile*  hold = TProfileHelper::RebinAxis(this, x, axis);
01512    if ( hold ) {
01513       fTsumwy  = hold->fTsumwy;
01514       fTsumwy2 = hold->fTsumwy2;
01515       
01516       delete hold;
01517    }
01518 }
01519 
01520 //______________________________________________________________________________
01521 void TProfile::Reset(Option_t *option)
01522 {
01523 //*-*-*-*-*-*-*-*-*-*Reset contents of a Profile histogram*-*-*-*-*-*-*-*-*
01524 //*-*                =====================================
01525    TH1D::Reset(option);
01526    fBinEntries.Reset();
01527    fBinSumw2.Reset();
01528    TString opt = option;
01529    opt.ToUpper();
01530    if (opt.Contains("ICE")) return;
01531    fTsumwy  = 0;
01532    fTsumwy2 = 0;
01533 }
01534 
01535 //______________________________________________________________________________
01536 void TProfile::SavePrimitive(ostream &out, Option_t *option /*= ""*/)
01537 {
01538     // Save primitive as a C++ statement(s) on output stream out
01539 
01540    //Note the following restrictions in the code generated:
01541    // - variable bin size not implemented
01542    // - SetErrorOption not implemented
01543 
01544    Bool_t nonEqiX = kFALSE;
01545    Int_t i;
01546    // Check if the profile has equidistant X bins or not.  If not, we
01547    // create an array holding the bins.
01548    if (GetXaxis()->GetXbins()->fN && GetXaxis()->GetXbins()->fArray) {
01549       nonEqiX = kTRUE;
01550       out << "   Double_t xAxis[" << GetXaxis()->GetXbins()->fN
01551           << "] = {";
01552       for (i = 0; i < GetXaxis()->GetXbins()->fN; i++) {
01553          if (i != 0) out << ", ";
01554          out << GetXaxis()->GetXbins()->fArray[i];
01555       }
01556       out << "}; " << endl;
01557    }
01558 
01559    char quote = '"';
01560    out<<"   "<<endl;
01561    out<<"   "<<ClassName()<<" *";
01562 
01563    //histogram pointer has by default teh histogram name.
01564    //however, in case histogram has no directory, it is safer to add a incremental suffix
01565    static Int_t hcounter = 0;
01566    TString histName = GetName();
01567    if (!fDirectory) {
01568       hcounter++;
01569       histName += "__";
01570       histName += hcounter;
01571    }
01572    const char *hname = histName.Data();
01573    
01574    out<<hname<<" = new "<<ClassName()<<"("<<quote<<GetName()<<quote<<","<<quote<<GetTitle()<<quote
01575                  <<","<<GetXaxis()->GetNbins();
01576    if (nonEqiX)
01577       out << ", xAxis";
01578    else
01579       out << "," << GetXaxis()->GetXmin()
01580           << "," << GetXaxis()->GetXmax()
01581           <<","<<quote<<GetErrorOption()<<quote<<");"<<endl;
01582 
01583    // save bin entries
01584    Int_t bin;
01585    for (bin=0;bin<fNcells;bin++) {
01586       Double_t bi = GetBinEntries(bin);
01587       if (bi) {
01588          out<<"   "<<hname<<"->SetBinEntries("<<bin<<","<<bi<<");"<<endl;
01589       }
01590    }
01591    //save bin contents
01592    for (bin=0;bin<fNcells;bin++) {
01593       Double_t bc = fArray[bin];
01594       if (bc) {
01595          out<<"   "<<hname<<"->SetBinContent("<<bin<<","<<bc<<");"<<endl;
01596       }
01597    }
01598    // save bin errors
01599    if (fSumw2.fN) {
01600       for (bin=0;bin<fNcells;bin++) {
01601          Double_t be = TMath::Sqrt(fSumw2.fArray[bin]);
01602          if (be) {
01603             out<<"   "<<hname<<"->SetBinError("<<bin<<","<<be<<");"<<endl;
01604          }
01605       }
01606    }
01607 
01608    TH1::SavePrimitiveHelp(out, hname, option);
01609 }
01610 
01611 //______________________________________________________________________________
01612 void TProfile::Scale(Double_t c1, Option_t * option)
01613 {
01614 // *-*-*-*-*Multiply this profile by a constant c1*-*-*-*-*-*-*-*-*
01615 // *-*      ======================================
01616 //
01617 //   this = c1*this
01618 //
01619 // This function uses the services of TProfile::Add
01620 //
01621    
01622    TProfileHelper::Scale(this, c1, option);
01623 }
01624 
01625 //______________________________________________________________________________
01626 void TProfile::SetBinEntries(Int_t bin, Double_t w)
01627 {
01628 //*-*-*-*-*-*-*-*-*Set the number of entries in bin*-*-*-*-*-*-*-*-*-*-*-*
01629 //*-*              ================================
01630 
01631    if (bin < 0 || bin >= fNcells) return;
01632    fBinEntries.fArray[bin] = w;
01633 }
01634 
01635 //______________________________________________________________________________
01636 void TProfile::SetBins(Int_t nx, Double_t xmin, Double_t xmax)
01637 {
01638 //*-*-*-*-*-*-*-*-*Redefine  x axis parameters*-*-*-*-*-*-*-*-*-*-*-*
01639 //*-*              ===========================
01640 
01641    fXaxis.Set(nx,xmin,xmax);
01642    fNcells = nx+2;
01643    SetBinsLength(fNcells);
01644    fBinEntries.Set(fNcells);
01645    fSumw2.Set(fNcells);
01646    if (fBinSumw2.fN) fBinSumw2.Set(fNcells);
01647 }
01648 
01649 //______________________________________________________________________________
01650 void TProfile::SetBins(Int_t nx, const Double_t *xbins)
01651 {
01652 //*-*-*-*-*-*-*-*-*Redefine  x axis parameters*-*-*-*-*-*-*-*-*-*-*-*
01653 //*-*              ===========================
01654 
01655    fXaxis.Set(nx,xbins);
01656    fNcells = nx+2;
01657    SetBinsLength(fNcells);
01658    fBinEntries.Set(fNcells);
01659    fSumw2.Set(fNcells);
01660    if (fBinSumw2.fN) fBinSumw2.Set(fNcells);
01661 }
01662 
01663 
01664 //______________________________________________________________________________
01665 void TProfile::SetBuffer(Int_t buffersize, Option_t *)
01666 {
01667 // set the buffer size in units of 8 bytes (double)
01668 
01669    if (fBuffer) {
01670       BufferEmpty();
01671       delete [] fBuffer;
01672       fBuffer = 0;
01673    }
01674    if (buffersize <= 0) {
01675       fBufferSize = 0;
01676       return;
01677    }
01678    if (buffersize < 100) buffersize = 100;
01679    fBufferSize = 1 + 3*buffersize;
01680    fBuffer = new Double_t[fBufferSize];
01681    memset(fBuffer,0,8*fBufferSize);
01682 }
01683 
01684 //______________________________________________________________________________
01685 void TProfile::SetErrorOption(Option_t *option)
01686 {
01687 //*-*-*-*-*-*-*-*-*-*Set option to compute profile errors*-*-*-*-*-*-*-*-*
01688 //*-*                =====================================
01689 //
01690 //    The computation of errors is based on the parameter option:
01691 //    option:
01692 //     ' '  (Default) Errors are Spread/SQRT(N) for Spread.ne.0. ,
01693 //                      "     "  SQRT(Y)/SQRT(N) for Spread.eq.0,N.gt.0 ,
01694 //                      "     "  0.  for N.eq.0
01695 //     's'            Errors are Spread  for Spread.ne.0. ,
01696 //                      "     "  SQRT(Y)  for Spread.eq.0,N.gt.0 ,
01697 //                      "     "  0.  for N.eq.0
01698 //     'i'            Errors are Spread/SQRT(N) for Spread.ne.0. ,
01699 //                      "     "  1./SQRT(12.*N) for Spread.eq.0,N.gt.0 ,
01700 //                      "     "  0.  for N.eq.0
01701 //     'g'            Errors are 1./SQRT(W) for Spread.ne.0. , 
01702 //                      "     "  0.  for N.eq.0
01703 //                    W is the sum of weights of the profile. 
01704 //                    This option is for measurements y +/ dy and  the profile is filled with 
01705 //                    weights w = 1/dy**2
01706 //
01707 //   See TProfile::BuildOptions for explanation of all options
01708 
01709    TString opt = option;
01710    opt.ToLower();
01711    fErrorMode = kERRORMEAN;
01712    if (opt.Contains("s")) fErrorMode = kERRORSPREAD;
01713    if (opt.Contains("i")) fErrorMode = kERRORSPREADI;
01714    if (opt.Contains("g")) fErrorMode = kERRORSPREADG;
01715 }
01716 
01717 //______________________________________________________________________________
01718 void TProfile::Streamer(TBuffer &R__b)
01719 {
01720    // Stream an object of class TProfile.
01721 
01722    if (R__b.IsReading()) {
01723       UInt_t R__s, R__c;
01724       Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
01725       if (R__v > 2) {
01726          R__b.ReadClassBuffer(TProfile::Class(), this, R__v, R__s, R__c);
01727          return;
01728       }
01729       //====process old versions before automatic schema evolution
01730       TH1D::Streamer(R__b);
01731       fBinEntries.Streamer(R__b);
01732       Int_t errorMode;
01733       R__b >> errorMode;
01734       fErrorMode = (EErrorType)errorMode;
01735       if (R__v < 2) {
01736          Float_t ymin,ymax;
01737          R__b >> ymin; fYmin = ymin;
01738          R__b >> ymax; fYmax = ymax;
01739       } else {
01740          R__b >> fYmin;
01741          R__b >> fYmax;
01742       }
01743       R__b.CheckByteCount(R__s, R__c, TProfile::IsA());
01744       //====end of old versions
01745 
01746    } else {
01747       R__b.WriteClassBuffer(TProfile::Class(),this);
01748    }
01749 }
01750 //______________________________________________________________________________
01751 void TProfile::Sumw2()
01752 {
01753    // Create structure to store sum of squares of weights per bin  *-*-*-*-*-*-*-*
01754    //   This is needed to compute  the correct statistical quantities  
01755    //    of a profile filled with weights 
01756    //  
01757    //
01758    //  This function is automatically called when the histogram is created
01759    //  if the static function TH1::SetDefaultSumw2 has been called before.
01760 
01761    TProfileHelper::Sumw2(this);
01762 }

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