TProfile3D.cxx

Go to the documentation of this file.
00001 // @(#)root/hist:$Id: TProfile3D.cxx 38022 2011-02-09 16:12:32Z moneta $
00002 // Author: Rene Brun   17/05/2006
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 "TProfile3D.h"
00013 #include "THashList.h"
00014 #include "TMath.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 TProfile3D::fgApproximate = kFALSE;
00024 
00025 ClassImp(TProfile3D)
00026 
00027 //______________________________________________________________________________
00028 //
00029 //  Profile3D histograms are used to display the mean
00030 //  value of T and its RMS for each cell in X,Y,Z.
00031 //  Profile3D histograms are in many cases an
00032 //  The inter-relation of three measured quantities X, Y, Z and T can always 
00033 //  be visualized by a four-dimensional histogram or scatter-plot; 
00034 //  its representation on the line-printer is not particularly
00035 //  satisfactory, except for sparse data. If T is an unknown (but single-valued)
00036 //  approximate function of X,Y,Z this function is displayed by a profile3D histogram with
00037 //  much better precision than by a scatter-plot.
00038 //
00039 //  The following formulae show the cumulated contents (capital letters) and the values
00040 //  displayed by the printing or plotting routines (small letters) of the elements for cell I, J.
00041 //
00042 //                                                                2
00043 //      H(I,J,K)  =  sum T                      E(I,J,K)  =  sum T
00044 //      l(I,J,K)  =  sum l                      L(I,J,K)  =  sum l
00045 //      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)
00046 //      e(I,J,K)  =  s(I,J,K)/sqrt(L(I,J,K))
00047 //
00048 //  In the special case where s(I,J,K) is zero (eg, case of 1 entry only in one cell)
00049 //  e(I,J,K) is computed from the average of the s(I,J,K) for all cells.
00050 //  This simple/crude approximation was suggested in order to keep the cell
00051 //  during a fit operation.
00052 //
00053 //           Example of a profile3D histogram
00054 //{
00055 //  TCanvas *c1 = new TCanvas("c1","Profile histogram example",200,10,700,500);
00056 //  hprof3d  = new TProfile3D("hprof3d","Profile of pt versus px, py and pz",40,-4,4,40,-4,4,40,0,20);
00057 //  Double_t px, py, pz, pt;
00058 //  TRandom3 r(0);
00059 //  for ( Int_t i=0; i<25000; i++) {
00060 //     r.Rannor(px,py);
00061 //     pz = px*px + py*py;
00062 //     pt = r.Landau(0,1);
00063 //     hprof3d->Fill(px,py,pz,pt,1);
00064 //  }
00065 //  hprof3d->Draw();
00066 //}
00067 //
00068 // NOTE: A TProfile3D is drawn as it was a simple TH3
00069 
00070 //______________________________________________________________________________
00071 TProfile3D::TProfile3D() : TH3D()
00072 {
00073 //*-*-*-*-*-*Default constructor for Profile3D histograms*-*-*-*-*-*-*-*-*
00074 //*-*        ============================================
00075    fTsumwt = fTsumwt2 = 0;
00076    fScaling = kFALSE;
00077    BuildOptions(0,0,"");
00078 }
00079 
00080 //______________________________________________________________________________
00081 TProfile3D::~TProfile3D()
00082 {
00083 //*-*-*-*-*-*Default destructor for Profile3D histograms*-*-*-*-*-*-*-*-*
00084 //*-*        ===========================================
00085 
00086 }
00087 
00088 //______________________________________________________________________________
00089 TProfile3D::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)
00090     : TH3D(name,title,nx,xlow,xup,ny,ylow,yup,nz,zlow,zup)
00091 {
00092 //*-*-*-*-*-*Normal Constructor for Profile histograms*-*-*-*-*-*-*-*-*-*
00093 //*-*        ==========================================
00094 //
00095 //  The first eleven parameters are similar to TH3D::TH3D.
00096 //  All values of t are accepted at filling time.
00097 //  To fill a profile3D histogram, one must use TProfile3D::Fill function.
00098 //
00099 //  Note that when filling the profile histogram the function Fill
00100 //  checks if the variable t is betyween fTmin and fTmax.
00101 //  If a minimum or maximum value is set for the T scale before filling,
00102 //  then all values below tmin or above tmax will be discarded.
00103 //  Setting the minimum or maximum value for the T scale before filling
00104 //  has the same effect as calling the special TProfile3D constructor below
00105 //  where tmin and tmax are specified.
00106 //
00107 //  H(I,J,K) is printed as the cell contents. The errors computed are s(I,J,K) if CHOPT='S'
00108 //  (spread option), or e(I,J,K) if CHOPT=' ' (error on mean).
00109 //
00110 //        See TProfile3D::BuildOptions for explanation of parameters
00111 //
00112 //   see other constructors below with all possible combinations of
00113 //   fix and variable bin size like in TH3D.
00114 
00115    BuildOptions(0,0,option);
00116    if (xlow >= xup || ylow >= yup || zlow >= zup) SetBuffer(fgBufferSize);
00117 }
00118 
00119 //______________________________________________________________________________
00120 TProfile3D::TProfile3D(const char *name,const char *title,Int_t nx,const Double_t *xbins,Int_t ny,const Double_t *ybins,Int_t nz,const Double_t *zbins,Option_t *option)
00121     : TH3D(name,title,nx,xbins,ny,ybins,nz,zbins)
00122 {
00123 //  Create a 3-D Profile with variable bins in X , Y and Z
00124 
00125    BuildOptions(0,0,option);
00126 }
00127 
00128 //______________________________________________________________________________
00129 void TProfile3D::BuildOptions(Double_t tmin, Double_t tmax, Option_t *option)
00130 {
00131 //*-*-*-*-*-*-*Set Profile3D histogram structure and options*-*-*-*-*-*-*-*-*
00132 //*-*          =============================================
00133 //
00134 //    If a cell has N data points all with the same value T (especially
00135 //    possible when dealing with integers), the spread in T for that cell
00136 //    is zero, and the uncertainty assigned is also zero, and the cell is
00137 //    ignored in making subsequent fits. If SQRT(T) was the correct error
00138 //    in the case above, then SQRT(T)/SQRT(N) would be the correct error here.
00139 //    In fact, any cell with non-zero number of entries N but with zero spread
00140 //    should have an uncertainty SQRT(T)/SQRT(N).
00141 //
00142 //    Now, is SQRT(T)/SQRT(N) really the correct uncertainty?
00143 //    that it is only in the case where the T variable is some sort
00144 //    of counting statistics, following a Poisson distribution. This should
00145 //    probably be set as the default case. However, T can be any variable
00146 //    from an original NTUPLE, not necessarily distributed "Poissonly".
00147 //    The computation of errors is based on the parameter option:
00148 //    option:
00149 //     ' '  (Default) Errors are Spread/SQRT(N) for Spread.ne.0. ,
00150 //                      "     "  SQRT(T)/SQRT(N) for Spread.eq.0,N.gt.0 ,
00151 //                      "     "  0.  for N.eq.0
00152 //     's'            Errors are Spread  for Spread.ne.0. ,
00153 //                      "     "  SQRT(T)  for Spread.eq.0,N.gt.0 ,
00154 //                      "     "  0.  for N.eq.0
00155 //     'i'            Errors are Spread/SQRT(N) for Spread.ne.0. ,
00156 //                      "     "  1./SQRT(12.*N) for Spread.eq.0,N.gt.0 ,
00157 //                      "     "  0.  for N.eq.0
00158 //
00159 //    The third case above corresponds to Integer T values for which the
00160 //    uncertainty is +-0.5, with the assumption that the probability that T
00161 //    takes any value between T-0.5 and T+0.5 is uniform (the same argument
00162 //    goes for T uniformly distributed between T and T+1); this would be
00163 //    useful if T is an ADC measurement, for example. Other, fancier options
00164 //    would be possible, at the cost of adding one more parameter to the PROFILE2D
00165 //    For example, if all T variables are distributed according to some
00166 //    known Gaussian of standard deviation Sigma, then:
00167 //     'G'            Errors are Spread/SQRT(N) for Spread.ne.0. ,
00168 //                      "     "  Sigma/SQRT(N) for Spread.eq.0,N.gt.0 ,
00169 //                      "     "  0.  for N.eq.0
00170 //    For example, this would be useful when all T's are experimental quantities
00171 //    measured with the same instrument with precision Sigma.
00172 //
00173 //
00174 
00175    SetErrorOption(option);
00176 
00177    fBinEntries.Set(fNcells);  //*-* create number of entries per cell array
00178 
00179    TH1::Sumw2();                   //*-* create sum of squares of weights array times y 
00180    if (fgDefaultSumw2) Sumw2();    // optionally create sum of squares of weights
00181 
00182    fTmin = tmin;
00183    fTmax = tmax;
00184    fScaling = kFALSE;
00185    fTsumwt  = fTsumwt2 = 0;
00186 }
00187 
00188 //______________________________________________________________________________
00189 TProfile3D::TProfile3D(const TProfile3D &profile) : TH3D()
00190 {
00191    //copy constructor
00192    ((TProfile3D&)profile).Copy(*this);
00193 }
00194 
00195 
00196 //______________________________________________________________________________
00197 void TProfile3D::Add(TF1 *, Double_t , Option_t*)
00198 {
00199    // Performs the operation: this = this + c1*f1
00200 
00201    Error("Add","Function not implemented for TProfile3D");
00202    return;
00203 }
00204 
00205 
00206 //______________________________________________________________________________
00207 void TProfile3D::Add(const TH1 *h1, Double_t c1)
00208 {
00209    // Performs the operation: this = this + c1*h1
00210 
00211    if (!h1) {
00212       Error("Add","Attempt to add a non-existing profile");
00213       return;
00214    }
00215    if (!h1->InheritsFrom(TProfile3D::Class())) {
00216       Error("Add","Attempt to add a non-profile2D object");
00217       return;
00218    }
00219 
00220    TProfileHelper::Add(this, this, h1, 1, c1);
00221 }
00222 
00223 //______________________________________________________________________________
00224 void TProfile3D::Add(const TH1 *h1, const TH1 *h2, Double_t c1, Double_t c2)
00225 {
00226 //*-*-*-*-*Replace contents of this profile3D by the addition of h1 and h2*-*-*
00227 //*-*      ===============================================================
00228 //
00229 //   this = c1*h1 + c2*h2
00230 //
00231 
00232    if (!h1 || !h2) {
00233       Error("Add","Attempt to add a non-existing profile");
00234       return;
00235    }
00236    if (!h1->InheritsFrom(TProfile3D::Class())) {
00237       Error("Add","Attempt to add a non-profile3D object");
00238       return;
00239    }
00240    if (!h2->InheritsFrom(TProfile3D::Class())) {
00241       Error("Add","Attempt to add a non-profile3D object");
00242       return;
00243    }
00244 
00245    TProfileHelper::Add(this, h1, h2, c1, c2);
00246 }
00247 
00248 
00249 //______________________________________________________________________________
00250 void TProfile3D::Approximate(Bool_t approx)
00251 {
00252 //     static function
00253 // set the fgApproximate flag. When the flag is true, the function GetBinError
00254 // will approximate the bin error with the average profile error on all bins
00255 // in the following situation only
00256 //  - the number of bins in the profile3D is less than 10404 (eg 100x100x100)
00257 //  - the bin number of entries is small ( <5)
00258 //  - the estimated bin error is extremely small compared to the bin content
00259 //  (see TProfile3D::GetBinError)
00260 
00261    fgApproximate = approx;
00262 }
00263 
00264 
00265 //______________________________________________________________________________
00266 Int_t TProfile3D::BufferEmpty(Int_t action)
00267 {
00268 // Fill histogram with all entries in the buffer.
00269 // action = -1 histogram is reset and refilled from the buffer (called by THistPainter::Paint)
00270 // action =  0 histogram is filled from the buffer
00271 // action =  1 histogram is filled and buffer is deleted
00272 //             The buffer is automatically deleted when the number of entries
00273 //             in the buffer is greater than the number of entries in the histogram
00274 
00275    // do we need to compute the bin size?
00276    if (!fBuffer) return 0;
00277    Int_t nbentries = (Int_t)fBuffer[0];
00278    if (!nbentries) return 0;
00279    Double_t *buffer = fBuffer;
00280    if (nbentries < 0) {
00281       if (action == 0) return 0;
00282       nbentries  = -nbentries;
00283       fBuffer=0;
00284       Reset();
00285       fBuffer = buffer;
00286    }
00287    if (TestBit(kCanRebin) || fXaxis.GetXmax() <= fXaxis.GetXmin() || fYaxis.GetXmax() <= fYaxis.GetXmin()) {
00288       //find min, max of entries in buffer
00289       Double_t xmin = fBuffer[2];
00290       Double_t xmax = xmin;
00291       Double_t ymin = fBuffer[3];
00292       Double_t ymax = ymin;
00293       Double_t zmin = fBuffer[4];
00294       Double_t zmax = zmin;
00295       for (Int_t i=1;i<nbentries;i++) {
00296          Double_t x = fBuffer[5*i+2];
00297          if (x < xmin) xmin = x;
00298          if (x > xmax) xmax = x;
00299          Double_t y = fBuffer[5*i+3];
00300          if (y < ymin) ymin = y;
00301          if (y > ymax) ymax = y;
00302          Double_t z = fBuffer[5*i+4];
00303          if (z < zmin) zmin = z;
00304          if (z > zmax) zmax = z;
00305      }
00306       if (fXaxis.GetXmax() <= fXaxis.GetXmin() || fYaxis.GetXmax() <= fYaxis.GetXmin() || fZaxis.GetXmax() <= fZaxis.GetXmin()) {
00307          THLimitsFinder::GetLimitsFinder()->FindGoodLimits(this,xmin,xmax,ymin,ymax,zmin,zmax);
00308       } else {
00309          fBuffer = 0;
00310          Int_t keep = fBufferSize; fBufferSize = 0;
00311          if (xmin <  fXaxis.GetXmin()) RebinAxis(xmin,&fXaxis);
00312          if (xmax >= fXaxis.GetXmax()) RebinAxis(xmax,&fXaxis);
00313          if (ymin <  fYaxis.GetXmin()) RebinAxis(ymin,&fYaxis);
00314          if (ymax >= fYaxis.GetXmax()) RebinAxis(ymax,&fYaxis);
00315          if (zmin <  fZaxis.GetXmin()) RebinAxis(zmin,&fZaxis);
00316          if (zmax >= fZaxis.GetXmax()) RebinAxis(zmax,&fZaxis);
00317          fBuffer = buffer;
00318          fBufferSize = keep;
00319       }
00320    }
00321 
00322    fBuffer = 0;
00323    for (Int_t i=0;i<nbentries;i++) {
00324       Fill(buffer[5*i+2],buffer[5*i+3],buffer[5*i+4],buffer[5*i+5],buffer[5*i+1]);
00325    }
00326    fBuffer = buffer;
00327 
00328    if (action > 0) { delete [] fBuffer; fBuffer = 0; fBufferSize = 0;}
00329    else {
00330       if (nbentries == (Int_t)fEntries) fBuffer[0] = -nbentries;
00331       else                              fBuffer[0] = 0;
00332    }
00333    return nbentries;
00334 }
00335 
00336 //______________________________________________________________________________
00337 Int_t TProfile3D::BufferFill(Double_t x, Double_t y, Double_t z, Double_t t, Double_t w)
00338 {
00339 // accumulate arguments in buffer. When buffer is full, empty the buffer
00340 // fBuffer[0] = number of entries in buffer
00341 // fBuffer[1] = w of first entry
00342 // fBuffer[2] = x of first entry
00343 // fBuffer[3] = y of first entry
00344 // fBuffer[4] = z of first entry
00345 // fBuffer[5] = t of first entry
00346 
00347    if (!fBuffer) return -3;
00348    Int_t nbentries = (Int_t)fBuffer[0];
00349    if (nbentries < 0) {
00350       nbentries  = -nbentries;
00351       fBuffer[0] =  nbentries;
00352       if (fEntries > 0) {
00353          Double_t *buffer = fBuffer; fBuffer=0;
00354          Reset();
00355          fBuffer = buffer;
00356       }
00357    }
00358    if (5*nbentries+5 >= fBufferSize) {
00359       BufferEmpty(1);
00360       return Fill(x,y,z,t,w);
00361    }
00362    fBuffer[5*nbentries+1] = w;
00363    fBuffer[5*nbentries+2] = x;
00364    fBuffer[5*nbentries+3] = y;
00365    fBuffer[5*nbentries+4] = z;
00366    fBuffer[5*nbentries+5] = t;
00367    fBuffer[0] += 1;
00368    return -2;
00369 }
00370 
00371 //______________________________________________________________________________
00372 void TProfile3D::Copy(TObject &obj) const
00373 {
00374 //*-*-*-*-*-*-*-*Copy a Profile3D histogram to a new profile2D histogram*-*-*-*
00375 //*-*            =======================================================
00376 
00377    TH3D::Copy(((TProfile3D&)obj));
00378    fBinEntries.Copy(((TProfile3D&)obj).fBinEntries);
00379    fBinSumw2.Copy(((TProfile3D&)obj).fBinSumw2);
00380    for (int bin=0;bin<fNcells;bin++) {
00381       ((TProfile3D&)obj).fArray[bin]        = fArray[bin];
00382       ((TProfile3D&)obj).fSumw2.fArray[bin] = fSumw2.fArray[bin];
00383    }
00384    ((TProfile3D&)obj).fTmin = fTmin;
00385    ((TProfile3D&)obj).fTmax = fTmax;
00386    ((TProfile3D&)obj).fScaling = fScaling;
00387    ((TProfile3D&)obj).fErrorMode = fErrorMode;
00388    ((TProfile3D&)obj).fTsumwt  = fTsumwt;
00389    ((TProfile3D&)obj).fTsumwt2 = fTsumwt2;
00390 }
00391 
00392 
00393 //______________________________________________________________________________
00394 void TProfile3D::Divide(TF1 *, Double_t )
00395 {
00396    // Performs the operation: this = this/(c1*f1)
00397 
00398    Error("Divide","Function not implemented for TProfile3D");
00399    return;
00400 }
00401 
00402 //______________________________________________________________________________
00403 void TProfile3D::Divide(const TH1 *h1)
00404 {
00405 //*-*-*-*-*-*-*-*-*-*-*Divide this profile2D by h1*-*-*-*-*-*-*-*-*-*-*-*-*
00406 //*-*                  ===========================
00407 //
00408 //   this = this/h1
00409 //
00410 
00411    if (!h1) {
00412       Error("Divide","Attempt to divide a non-existing profile2D");
00413       return;
00414    }
00415    if (!h1->InheritsFrom(TProfile3D::Class())) {
00416       Error("Divide","Attempt to divide a non-profile3D object");
00417       return;
00418    }
00419    TProfile3D *p1 = (TProfile3D*)h1;
00420 
00421 //*-*- Check profile compatibility
00422    Int_t nx = GetNbinsX();
00423    if (nx != p1->GetNbinsX()) {
00424       Error("Divide","Attempt to divide profiles with different number of bins");
00425       return;
00426    }
00427    Int_t ny = GetNbinsY();
00428    if (ny != p1->GetNbinsY()) {
00429       Error("Divide","Attempt to divide profiles with different number of bins");
00430       return;
00431    }
00432    Int_t nz = GetNbinsZ();
00433    if (nz != p1->GetNbinsZ()) {
00434       Error("Divide","Attempt to divide profiles with different number of bins");
00435       return;
00436    }
00437 
00438 //*-*- Reset statistics
00439    fEntries = fTsumw   = fTsumw2 = fTsumwx = fTsumwx2 = 0;
00440 
00441 //*-*- Loop on bins (including underflows/overflows)
00442    Int_t bin,binx,biny,binz;
00443    Double_t *cu1 = p1->GetW();
00444    Double_t *er1 = p1->GetW2();
00445    Double_t *en1 = p1->GetB();
00446    Double_t c0,c1,w,u,x,y,z;
00447    for (binx =0;binx<=nx+1;binx++) {
00448       for (biny =0;biny<=ny+1;biny++) {
00449          for (binz =0;binz<=nz+1;binz++) {
00450             bin   = GetBin(binx,biny,binz);
00451             c0  = fArray[bin];
00452             c1  = cu1[bin];
00453             if (c1) w = c0/c1;
00454             else    w = 0;
00455             fArray[bin] = w;
00456             u = TMath::Abs(w);
00457             x = fXaxis.GetBinCenter(binx);
00458             y = fYaxis.GetBinCenter(biny);
00459             z = fZaxis.GetBinCenter(binz);
00460             fEntries++;
00461             fTsumw   += u;
00462             fTsumw2  += u*u;
00463             fTsumwx  += u*x;
00464             fTsumwx2 += u*x*x;
00465             fTsumwy  += u*y;
00466             fTsumwy2 += u*y*y;
00467             fTsumwxy += u*x*y;
00468             fTsumwz  += u;
00469             fTsumwz2 += u*z;
00470             fTsumwxz += u*x*z;
00471             fTsumwyz += u*y*z;
00472             fTsumwt  += u;
00473             fTsumwt2 += u*u;
00474             Double_t e0 = fSumw2.fArray[bin];
00475             Double_t e1 = er1[bin];
00476             Double_t c12= c1*c1;
00477             if (!c1) fSumw2.fArray[bin] = 0;
00478             else     fSumw2.fArray[bin] = (e0*c1*c1 + e1*c0*c0)/(c12*c12);
00479             if (!en1[bin]) fBinEntries.fArray[bin] = 0;
00480             else           fBinEntries.fArray[bin] /= en1[bin];
00481          }
00482       }
00483    }
00484    // mantaining the correct sum of weights square is not supported when dividing
00485    // bin error resulting from division of profile needs to be checked 
00486    if (fBinSumw2.fN) { 
00487       Warning("Divide","Cannot preserve during the division of profiles the sum of bin weight square");
00488       fBinSumw2 = TArrayD();
00489    }
00490 }
00491 
00492 
00493 //______________________________________________________________________________
00494 void TProfile3D::Divide(const TH1 *h1, const TH1 *h2, Double_t c1, Double_t c2, Option_t *option)
00495 {
00496 //*-*-*-*-*Replace contents of this profile2D by the division of h1 by h2*-*-*
00497 //*-*      ==============================================================
00498 //
00499 //   this = c1*h1/(c2*h2)
00500 //
00501 
00502    TString opt = option;
00503    opt.ToLower();
00504    Bool_t binomial = kFALSE;
00505    if (opt.Contains("b")) binomial = kTRUE;
00506    if (!h1 || !h2) {
00507       Error("Divide","Attempt to divide a non-existing profile2D");
00508       return;
00509    }
00510    if (!h1->InheritsFrom(TProfile3D::Class())) {
00511       Error("Divide","Attempt to divide a non-profile2D object");
00512       return;
00513    }
00514    TProfile3D *p1 = (TProfile3D*)h1;
00515    if (!h2->InheritsFrom(TProfile3D::Class())) {
00516       Error("Divide","Attempt to divide a non-profile2D object");
00517       return;
00518    }
00519    TProfile3D *p2 = (TProfile3D*)h2;
00520 
00521 //*-*- Check histogram compatibility
00522    Int_t nx = GetNbinsX();
00523    if (nx != p1->GetNbinsX() || nx != p2->GetNbinsX()) {
00524       Error("Divide","Attempt to divide profiles with different number of bins");
00525       return;
00526    }
00527    Int_t ny = GetNbinsY();
00528    if (ny != p1->GetNbinsY() || ny != p2->GetNbinsY()) {
00529       Error("Divide","Attempt to divide profiles with different number of bins");
00530       return;
00531    }
00532    Int_t nz = GetNbinsZ();
00533    if (nz != p1->GetNbinsZ() || nz != p2->GetNbinsZ()) {
00534       Error("Divide","Attempt to divide profiles with different number of bins");
00535       return;
00536    }
00537    if (!c2) {
00538       Error("Divide","Coefficient of dividing profile cannot be zero");
00539       return;
00540    }
00541 
00542 //*-*- Reset statistics
00543    fEntries = fTsumw   = fTsumw2 = fTsumwx = fTsumwx2 = 0;
00544 
00545 //*-*- Loop on bins (including underflows/overflows)
00546    Int_t bin,binx,biny,binz;
00547    Double_t *cu1 = p1->GetW();
00548    Double_t *cu2 = p2->GetW();
00549    Double_t *er1 = p1->GetW2();
00550    Double_t *er2 = p2->GetW2();
00551    Double_t *en1 = p1->GetB();
00552    Double_t *en2 = p2->GetB();
00553    Double_t b1,b2,w,u,x,y,z,ac1,ac2;
00554    ac1 = TMath::Abs(c1);
00555    ac2 = TMath::Abs(c2);
00556    for (binx =0;binx<=nx+1;binx++) {
00557       for (biny =0;biny<=ny+1;biny++) {
00558          for (binz =0;binz<=nz+1;binz++) {
00559             bin   = GetBin(binx,biny,binz);
00560             b1  = cu1[bin];
00561             b2  = cu2[bin];
00562             if (b2) w = c1*b1/(c2*b2);
00563             else    w = 0;
00564             fArray[bin] = w;
00565             u = TMath::Abs(w);
00566             x = fXaxis.GetBinCenter(binx);
00567             y = fYaxis.GetBinCenter(biny);
00568             z = fZaxis.GetBinCenter(biny);
00569             fEntries++;
00570             fTsumw   += u;
00571             fTsumw2  += u*u;
00572             fTsumwx  += u*x;
00573             fTsumwx2 += u*x*x;
00574             fTsumwy  += u*y;
00575             fTsumwy2 += u*y*y;
00576             fTsumwxy += u*x*y;
00577             fTsumwz  += u*z;
00578             fTsumwz2 += u*z*z;
00579             fTsumwxz += u*x*z;
00580             fTsumwyz += u*y*z;
00581             fTsumwt  += u;
00582             fTsumwt2 += u*u;
00583             Double_t e1 = er1[bin];
00584             Double_t e2 = er2[bin];
00585           //Double_t b22= b2*b2*d2;
00586             Double_t b22= b2*b2*TMath::Abs(c2);
00587             if (!b2) fSumw2.fArray[bin] = 0;
00588             else {
00589                if (binomial) {
00590                   fSumw2.fArray[bin] = TMath::Abs(w*(1-w)/(c2*b2));
00591                } else {
00592                   fSumw2.fArray[bin] = ac1*ac2*(e1*b2*b2 + e2*b1*b1)/(b22*b22);
00593                }
00594             }
00595             if (!en2[bin]) fBinEntries.fArray[bin] = 0;
00596             else           fBinEntries.fArray[bin] = en1[bin]/en2[bin];
00597          }
00598       }
00599    }
00600 }
00601 
00602 //______________________________________________________________________________
00603 TH1 *TProfile3D::DrawCopy(Option_t *option) const
00604 {
00605 //*-*-*-*-*-*-*-*Draw a copy of this profile3D histogram*-*-*-*-*-*-*-*-*-*-*
00606 //*-*            =======================================
00607    TString opt = option;
00608    opt.ToLower();
00609    if (gPad && !opt.Contains("same")) gPad->Clear();
00610    TProfile3D *newpf = new TProfile3D();
00611    Copy(*newpf);
00612    newpf->SetDirectory(0);
00613    newpf->SetBit(kCanDelete);
00614    newpf->AppendPad(option);
00615    return newpf;
00616 }
00617 
00618 //______________________________________________________________________________
00619 Int_t TProfile3D::Fill(Double_t x, Double_t y, Double_t z, Double_t t)
00620 {
00621 //*-*-*-*-*-*-*-*-*-*-*Fill a Profile3D histogram (no weights)*-*-*-*-*-*-*-*
00622 //*-*                  =======================================
00623 
00624    if (fBuffer) return BufferFill(x,y,z,t,1);
00625 
00626    Int_t bin,binx,biny,binz;
00627 
00628    if (fTmin != fTmax) {
00629       if (t <fTmin || t> fTmax) return -1;
00630    }
00631 
00632    fEntries++;
00633    binx =fXaxis.FindBin(x);
00634    biny =fYaxis.FindBin(y);
00635    binz =fZaxis.FindBin(z);
00636    if (binx <0 || biny <0 || binz<0) return -1;
00637    bin  = GetBin(binx,biny,binz);
00638    AddBinContent(bin, t);
00639    fSumw2.fArray[bin] += (Double_t)t*t;
00640    fBinEntries.fArray[bin] += 1;
00641    if (fBinSumw2.fN)  fBinSumw2.fArray[bin] += 1;
00642    if (binx == 0 || binx > fXaxis.GetNbins()) {
00643       if (!fgStatOverflows) return -1;
00644    }
00645    if (biny == 0 || biny > fYaxis.GetNbins()) {
00646       if (!fgStatOverflows) return -1;
00647    }
00648    if (binz == 0 || binz > fZaxis.GetNbins()) {
00649       if (!fgStatOverflows) return -1;
00650    }
00651 //printf("x=%g, y=%g, z=%g, t=%g, binx=%d, biny=%d, binz=%d, bin=%d\n",x,y,z,t,binx,biny,binz,bin);
00652    ++fTsumw;
00653    ++fTsumw2;
00654    fTsumwx  += x;
00655    fTsumwx2 += x*x;
00656    fTsumwy  += y;
00657    fTsumwy2 += y*y;
00658    fTsumwxy += x*y;
00659    fTsumwz  += z;
00660    fTsumwz2 += z*z;
00661    fTsumwxz += x*z;
00662    fTsumwyz += y*z;
00663    fTsumwt  += t;
00664    fTsumwt2 += t*t;
00665    return bin;
00666 }
00667 
00668 //______________________________________________________________________________
00669 Int_t TProfile3D::Fill(Double_t x, Double_t y, Double_t z, Double_t t, Double_t w)
00670 {
00671 //*-*-*-*-*-*-*-*-*-*-*Fill a Profile3D histogram with weights*-*-*-*-*-*-*-*
00672 //*-*                  =======================================
00673 
00674    if (fBuffer) return BufferFill(x,y,z,t,w);
00675 
00676    Int_t bin,binx,biny,binz;
00677 
00678    if (fTmin != fTmax) {
00679       if (t <fTmin || z> fTmax) return -1;
00680    }
00681 
00682    Double_t u= (w > 0 ? w : -w);
00683    fEntries++;
00684    binx =fXaxis.FindBin(x);
00685    biny =fYaxis.FindBin(y);
00686    binz =fZaxis.FindBin(z);
00687    if (binx <0 || biny <0 || binz<0) return -1;
00688    bin  = GetBin(binx,biny,binz);
00689    AddBinContent(bin, u*t);
00690    fSumw2.fArray[bin] += u*t*t;
00691    fBinEntries.fArray[bin] += u;
00692    if (fBinSumw2.fN)  fBinSumw2.fArray[bin] += u*u;
00693    if (binx == 0 || binx > fXaxis.GetNbins()) {
00694       if (!fgStatOverflows) return -1;
00695    }
00696    if (biny == 0 || biny > fYaxis.GetNbins()) {
00697       if (!fgStatOverflows) return -1;
00698    }
00699    if (binz == 0 || binz > fZaxis.GetNbins()) {
00700       if (!fgStatOverflows) return -1;
00701    }
00702    fTsumw   += u;
00703    fTsumw2  += u*u;
00704    fTsumwx  += u*x;
00705    fTsumwx2 += u*x*x;
00706    fTsumwy  += u*y;
00707    fTsumwy2 += u*y*y;
00708    fTsumwxy += u*x*y;
00709    fTsumwz  += u*z;
00710    fTsumwz2 += u*z*z;
00711    fTsumwxz += u*x*z;
00712    fTsumwyz += u*y*z;
00713    fTsumwt  += u*t;
00714    fTsumwt2 += u*t*t;
00715    return bin;
00716 }
00717 
00718 //______________________________________________________________________________
00719 Double_t TProfile3D::GetBinContent(Int_t bin) const
00720 {
00721 //*-*-*-*-*-*-*Return bin content of a Profile3D histogram*-*-*-*-*-*-*-*-*
00722 //*-*          ===========================================
00723 
00724    if (fBuffer) ((TProfile3D*)this)->BufferEmpty();
00725 
00726    if (bin < 0 || bin >= fNcells) return 0;
00727    if (fBinEntries.fArray[bin] == 0) return 0;
00728    if (!fArray) return 0;
00729    return fArray[bin]/fBinEntries.fArray[bin];
00730 }
00731 
00732 //______________________________________________________________________________
00733 Double_t TProfile3D::GetBinEntries(Int_t bin) const
00734 {
00735 //*-*-*-*-*-*-*Return bin entries of a Profile3D histogram*-*-*-*-*-*-*-*-*
00736 //*-*          ===========================================
00737 
00738    if (fBuffer) ((TProfile3D*)this)->BufferEmpty();
00739 
00740    if (bin < 0 || bin >= fNcells) return 0;
00741    return fBinEntries.fArray[bin];
00742 }
00743 
00744 //______________________________________________________________________________
00745 Double_t TProfile3D::GetBinEffectiveEntries(Int_t bin)
00746 {
00747 //            Return bin effective entries for a weighted filled Profile histogram. 
00748 //            In case of an unweighted profile, it is equivalent to the number of entries per bin   
00749 //            The effective entries is defined as the square of the sum of the weights divided by the 
00750 //            sum of the weights square. 
00751 //            TProfile::Sumw2() must be called before filling the profile with weights. 
00752 //            Only by calling this method the  sum of the square of the weights per bin is stored. 
00753 //  
00754 //*-*          =========================================
00755 
00756    return TProfileHelper::GetBinEffectiveEntries((TProfile3D*)this, bin);
00757 }
00758 
00759 //______________________________________________________________________________
00760 Double_t TProfile3D::GetBinError(Int_t bin) const
00761 {
00762 // *-*-*-*-*-*-*Return bin error of a Profile3D histogram*-*-*-*-*-*-*-*-*
00763 //
00764 // Computing errors: A moving field
00765 // =================================
00766 // The computation of errors for a TProfile3D has evolved with the versions
00767 // of ROOT. The difficulty is in computing errors for bins with low statistics.
00768 // - prior to version 3.10, we had no special treatment of low statistic bins.
00769 //   As a result, these bins had huge errors. The reason is that the
00770 //   expression eprim2 is very close to 0 (rounding problems) or 0.
00771 // - The algorithm is modified/protected for the case
00772 //   when a TProfile3D is projected (ProjectionX). The previous algorithm
00773 //   generated a N^2 problem when projecting a TProfile3D with a large number of
00774 //   bins (eg 100000).
00775 // - in version 3.10/02, a new static function TProfile::Approximate
00776 //   is introduced to enable or disable (default) the approximation.
00777 //   (see also comments in TProfile::GetBinError)
00778 
00779    return TProfileHelper::GetBinError((TProfile3D*)this, bin);
00780 }
00781 
00782 //______________________________________________________________________________
00783 Option_t *TProfile3D::GetErrorOption() const
00784 {
00785 //*-*-*-*-*-*-*-*-*-*Return option to compute profile2D errors*-*-*-*-*-*-*-*
00786 //*-*                =========================================
00787 
00788    if (fErrorMode == kERRORSPREAD)  return "s";
00789    if (fErrorMode == kERRORSPREADI) return "i";
00790    if (fErrorMode == kERRORSPREADG) return "g";
00791    return "";
00792 }
00793 
00794 //______________________________________________________________________________
00795 void TProfile3D::GetStats(Double_t *stats) const
00796 {
00797    // fill the array stats from the contents of this profile
00798    // The array stats must be correctly dimensionned in the calling program.
00799    // stats[0] = sumw
00800    // stats[1] = sumw2
00801    // stats[2] = sumwx
00802    // stats[3] = sumwx2
00803    // stats[4] = sumwy
00804    // stats[5] = sumwy2
00805    // stats[6] = sumwxy
00806    // stats[7] = sumwz
00807    // stats[8] = sumwz2
00808    // stats[9] = sumwxz
00809    // stats[10]= sumwyz
00810    // stats[11]= sumwt
00811    // stats[12]= sumwt2
00812    //
00813    // If no axis-subrange is specified (via TAxis::SetRange), the array stats
00814    // is simply a copy of the statistics quantities computed at filling time.
00815    // If a sub-range is specified, the function recomputes these quantities
00816    // from the bin contents in the current axis range.
00817 
00818    if (fBuffer) ((TProfile3D*)this)->BufferEmpty();
00819 
00820    // Loop on bins
00821    if (fTsumw == 0 || fXaxis.TestBit(TAxis::kAxisRange) || fYaxis.TestBit(TAxis::kAxisRange)) {
00822       Int_t bin, binx, biny,binz;
00823       Double_t w, w2;
00824       Double_t x,y,z;
00825       for (bin=0;bin<kNstat;bin++) stats[bin] = 0;
00826       if (!fBinEntries.fArray) return;
00827       for (binz=fZaxis.GetFirst();binz<=fZaxis.GetLast();binz++) {
00828          z = fZaxis.GetBinCenter(binz);
00829          for (biny=fYaxis.GetFirst();biny<=fYaxis.GetLast();biny++) {
00830             y = fYaxis.GetBinCenter(biny);
00831             for (binx=fXaxis.GetFirst();binx<=fXaxis.GetLast();binx++) {
00832                bin = GetBin(binx,biny,binz);
00833                w         = fBinEntries.fArray[bin];
00834                w2        = (fBinSumw2.fN ? fBinSumw2.fArray[bin] : w );
00835                x         = fXaxis.GetBinCenter(binx);
00836                stats[0]  += w;
00837                stats[1]  += w2;
00838                stats[2]  += w*x;
00839                stats[3]  += w*x*x;
00840                stats[4]  += w*y;
00841                stats[5]  += w*y*y;
00842                stats[6]  += w*x*y;
00843                stats[7]  += w*z;
00844                stats[8]  += w*z*z;
00845                stats[9]  += w*x*z;
00846                stats[10] += w*y*z;
00847                stats[11] += fArray[bin];
00848                stats[12] += fSumw2.fArray[bin];
00849             }
00850          }
00851       }
00852    } else {
00853       stats[0]  = fTsumw;
00854       stats[1]  = fTsumw2;
00855       stats[2]  = fTsumwx;
00856       stats[3]  = fTsumwx2;
00857       stats[4]  = fTsumwy;
00858       stats[5]  = fTsumwy2;
00859       stats[6]  = fTsumwxy;
00860       stats[7]  = fTsumwz;
00861       stats[8]  = fTsumwz2;
00862       stats[9]  = fTsumwxz;
00863       stats[10] = fTsumwyz;
00864       stats[11] = fTsumwt;
00865       stats[12] = fTsumwt2;
00866    }
00867 }
00868 
00869 //______________________________________________________________________________
00870 Long64_t TProfile3D::Merge(TCollection *li)
00871 {
00872    //Merge all histograms in the collection in this histogram.
00873    //This function computes the min/max for the axes,
00874    //compute a new number of bins, if necessary,
00875    //add bin contents, errors and statistics.
00876    //If overflows are present and limits are different the function will fail.
00877    //The function returns the total number of entries in the result histogram
00878    //if the merge is successfull, -1 otherwise.
00879    //
00880    //IMPORTANT remark. The 2 axis x and y may have different number
00881    //of bins and different limits, BUT the largest bin width must be
00882    //a multiple of the smallest bin width and the upper limit must also
00883    //be a multiple of the bin width.
00884 
00885    return TProfileHelper::Merge(this, li);
00886 
00887 //    if (!li) return 0;
00888 //    if (li->IsEmpty()) return (Int_t) GetEntries();
00889 
00890 //    TList inlist;
00891 //    TH1* hclone = (TH1*)Clone("FirstClone");
00892 //    R__ASSERT(hclone);
00893 //    BufferEmpty(1);         // To remove buffer.
00894 //    Reset();                // BufferEmpty sets limits so we can't use it later.
00895 //    SetEntries(0);
00896 //    inlist.Add(hclone);
00897 //    inlist.AddAll(li);
00898 
00899 //    TAxis newXAxis;
00900 //    TAxis newYAxis;
00901 //    TAxis newZAxis;
00902 //    Bool_t initialLimitsFound = kFALSE;
00903 //    Bool_t same = kTRUE;
00904 //    Bool_t allHaveLimits = kTRUE;
00905 
00906 //    TIter next(&inlist);
00907 //    while (TObject *o = next()) {
00908 //       TProfile3D* h = dynamic_cast<TProfile3D*> (o);
00909 //       if (!h) {
00910 //          Error("Add","Attempt to add object of class: %s to a %s",
00911 //                o->ClassName(),this->ClassName());
00912 //          return -1;
00913 //       }
00914 //       Bool_t hasLimits = h->GetXaxis()->GetXmin() < h->GetXaxis()->GetXmax();
00915 //       allHaveLimits = allHaveLimits && hasLimits;
00916 
00917 //       if (hasLimits) {
00918 //          h->BufferEmpty();
00919 //          if (!initialLimitsFound) {
00920 //             initialLimitsFound = kTRUE;
00921 //             newXAxis.Set(h->GetXaxis()->GetNbins(), h->GetXaxis()->GetXmin(),
00922 //                      h->GetXaxis()->GetXmax());
00923 //             newYAxis.Set(h->GetYaxis()->GetNbins(), h->GetYaxis()->GetXmin(),
00924 //                      h->GetYaxis()->GetXmax());
00925 //             newZAxis.Set(h->GetZaxis()->GetNbins(), h->GetZaxis()->GetXmin(),
00926 //                      h->GetZaxis()->GetXmax());
00927 //          }
00928 //          else {
00929 //             if (!RecomputeAxisLimits(newXAxis, *(h->GetXaxis()))) {
00930 //                Error("Merge", "Cannot merge histograms - limits are inconsistent:\n "
00931 //                      "first: (%d, %f, %f), second: (%d, %f, %f)",
00932 //                      newXAxis.GetNbins(), newXAxis.GetXmin(), newXAxis.GetXmax(),
00933 //                      h->GetXaxis()->GetNbins(), h->GetXaxis()->GetXmin(),
00934 //                      h->GetXaxis()->GetXmax());
00935 //             }
00936 //             if (!RecomputeAxisLimits(newYAxis, *(h->GetYaxis()))) {
00937 //                Error("Merge", "Cannot merge histograms - limits are inconsistent:\n "
00938 //                      "first: (%d, %f, %f), second: (%d, %f, %f)",
00939 //                      newYAxis.GetNbins(), newYAxis.GetXmin(), newYAxis.GetXmax(),
00940 //                      h->GetYaxis()->GetNbins(), h->GetYaxis()->GetXmin(),
00941 //                      h->GetYaxis()->GetXmax());
00942 //             }
00943 //             if (!RecomputeAxisLimits(newZAxis, *(h->GetZaxis()))) {
00944 //                Error("Merge", "Cannot merge histograms - limits are inconsistent:\n "
00945 //                      "first: (%d, %f, %f), second: (%d, %f, %f)",
00946 //                      newZAxis.GetNbins(), newZAxis.GetXmin(), newZAxis.GetXmax(),
00947 //                      h->GetZaxis()->GetNbins(), h->GetZaxis()->GetXmin(),
00948 //                      h->GetZaxis()->GetXmax());
00949 //             }
00950 //          }
00951 //       }
00952 //    }
00953 //    next.Reset();
00954 
00955 //    same = same && SameLimitsAndNBins(newXAxis, *GetXaxis())
00956 //                && SameLimitsAndNBins(newYAxis, *GetYaxis())
00957 //                && SameLimitsAndNBins(newZAxis, *GetZaxis());
00958 //    if (!same && initialLimitsFound)
00959 //       SetBins(newXAxis.GetNbins(), newXAxis.GetXmin(), newXAxis.GetXmax(),
00960 //               newYAxis.GetNbins(), newYAxis.GetXmin(), newYAxis.GetXmax(),
00961 //               newZAxis.GetNbins(), newZAxis.GetXmin(), newZAxis.GetXmax());
00962 
00963 //    if (!allHaveLimits) {
00964 //       // fill this histogram with all the data from buffers of histograms without limits
00965 //       while (TProfile3D* h = dynamic_cast<TProfile3D*> (next())) {
00966 //          if (h->GetXaxis()->GetXmin() >= h->GetXaxis()->GetXmax() && h->fBuffer) {
00967 //              // no limits
00968 //             Int_t nbentries = (Int_t)h->fBuffer[0];
00969 //             for (Int_t i = 0; i < nbentries; i++)
00970 //                Fill(h->fBuffer[5*i + 2], h->fBuffer[5*i + 3],
00971 //                     h->fBuffer[5*i + 4], h->fBuffer[5*i + 5], h->fBuffer[5*i + 1]);
00972 //          }
00973 //       }
00974 //       if (!initialLimitsFound)
00975 //          return (Int_t) GetEntries();  // all histograms have been processed
00976 //       next.Reset();
00977 //    }
00978 
00979 //    //merge bin contents and errors
00980 //    Double_t stats[kNstat], totstats[kNstat];
00981 //    for (Int_t i=0;i<kNstat;i++) {totstats[i] = stats[i] = 0;}
00982 //    GetStats(totstats);
00983 //    Double_t nentries = GetEntries();
00984 //    Int_t binx, biny, binz, ix, iy, iz, nx, ny, nz, bin, ibin;
00985 //    Int_t nbix = fXaxis.GetNbins();
00986 //    Int_t nbiy = fYaxis.GetNbins();
00987 //    Bool_t canRebin=TestBit(kCanRebin);
00988 //    ResetBit(kCanRebin); // reset, otherwise setting the under/overflow will rebin
00989 
00990 //    while (TProfile3D* h=(TProfile3D*)next()) {
00991 //       // process only if the histogram has limits; otherwise it was processed before
00992 //       if (h->GetXaxis()->GetXmin() < h->GetXaxis()->GetXmax()) {
00993 //          // import statistics
00994 //          h->GetStats(stats);
00995 //          for (Int_t i = 0; i < kNstat; i++)
00996 //             totstats[i] += stats[i];
00997 //          nentries += h->GetEntries();
00998 
00999 //          nx = h->GetXaxis()->GetNbins();
01000 //          ny = h->GetYaxis()->GetNbins();
01001 //          nz = h->GetZaxis()->GetNbins();
01002 
01003 //          for (binz = 0; binz <= nz + 1; binz++) {
01004 //             iz = fZaxis.FindBin(h->GetZaxis()->GetBinCenter(binz));
01005 //             for (biny = 0; biny <= ny + 1; biny++) {
01006 //                iy = fYaxis.FindBin(h->GetYaxis()->GetBinCenter(biny));
01007 //                for (binx = 0; binx <= nx + 1; binx++) {
01008 //                   ix = fXaxis.FindBin(h->GetXaxis()->GetBinCenter(binx));
01009 //                   bin  = binx +(nx+2)*(biny + (ny+2)*binz);
01010 //                   ibin = ix   +(nbix+2)*(iy + (nbiy+2)*iz);
01011 //                   if ((!same) && (binx == 0 || binx == nx + 1
01012 //                                || biny == 0 || biny == ny + 1 
01013 //                                || binz == 0 || binz == nz + 1)) {
01014 //                      if (h->GetW()[bin] != 0) {
01015 //                         Error("Merge", "Cannot merge histograms - the histograms have"
01016 //                                        " different limits and undeflows/overflows are present."
01017 //                                        " The initial histogram is now broken!");
01018 //                         return -1;
01019 //                      }
01020 //                   }
01021 //                   fArray[ibin]             += h->GetW()[bin];
01022 //                   fSumw2.fArray[ibin]      += h->GetW2()[bin];
01023 //                   fBinEntries.fArray[ibin] += h->GetB()[bin];
01024 //                }
01025 //             }
01026 //          }
01027 //          fEntries += h->GetEntries();
01028 //          fTsumw   += h->fTsumw;
01029 //          fTsumw2  += h->fTsumw2;
01030 //          fTsumwx  += h->fTsumwx;
01031 //          fTsumwx2 += h->fTsumwx2;
01032 //          fTsumwy  += h->fTsumwy;
01033 //          fTsumwy2 += h->fTsumwy2;
01034 //          fTsumwxy += h->fTsumwxy;
01035 //          fTsumwz  += h->fTsumwz;
01036 //          fTsumwz2 += h->fTsumwz2;
01037 //          fTsumwxz += h->fTsumwxz;
01038 //          fTsumwyz += h->fTsumwyz;
01039 //          fTsumwt  += h->fTsumwt;
01040 //          fTsumwt2 += h->fTsumwt2;
01041 //       }
01042 //    }
01043 //    if (canRebin) SetBit(kCanRebin);
01044 
01045 //    //copy merged stats
01046 //    PutStats(totstats);
01047 //    SetEntries(nentries);
01048 //    inlist.Remove(hclone);
01049 //    delete hclone;
01050 //    return (Long64_t)nentries;
01051 }
01052 
01053 //______________________________________________________________________________
01054 void TProfile3D::Multiply(TF1 *, Double_t )
01055 {
01056    // Performs the operation: this = this*c1*f1
01057 
01058    Error("Multiply","Function not implemented for TProfile3D");
01059    return;
01060 }
01061 
01062 //______________________________________________________________________________
01063 void TProfile3D::Multiply(const TH1 *)
01064 {
01065 //*-*-*-*-*-*-*-*-*-*-*Multiply this profile2D by h1*-*-*-*-*-*-*-*-*-*-*-*
01066 //*-*                  =============================
01067 //
01068 //   this = this*h1
01069 //
01070    Error("Multiply","Multiplication of profile2D histograms not implemented");
01071 }
01072 
01073 
01074 //______________________________________________________________________________
01075 void TProfile3D::Multiply(const TH1 *, const TH1 *, Double_t, Double_t, Option_t *)
01076 {
01077 //*-*-*-*-*Replace contents of this profile2D by multiplication of h1 by h2*-*
01078 //*-*      ================================================================
01079 //
01080 //   this = (c1*h1)*(c2*h2)
01081 //
01082 
01083    Error("Multiply","Multiplication of profile2D histograms not implemented");
01084 }
01085 
01086 //______________________________________________________________________________
01087 TH3D *TProfile3D::ProjectionXYZ(const char *name, Option_t *option) const
01088 {
01089 //*-*-*-*-*Project this profile3D into a 3-D histogram along X,Y,Z*-*-*-*-*-*-*
01090 //*-*      =====================================================
01091 //
01092 //   The projection is always of the type TH3D.
01093 //
01094 //   if option "E" is specified, the errors are computed. (default)
01095 //   if option "B" is specified, the content of bin of the returned histogram
01096 //      will be equal to the GetBinEntries(bin) of the profile,
01097 //   if option "C=E" the bin contents of the projection are set to the
01098 //       bin errors of the profile
01099 //
01100 
01101    TString opt = option;
01102    opt.ToLower();
01103    Int_t nx = fXaxis.GetNbins();
01104    Int_t ny = fYaxis.GetNbins();
01105    Int_t nz = fZaxis.GetNbins();
01106 
01107    // Create the projection histogram
01108    TString pname = name; 
01109    if (pname == "_px") { 
01110       pname = GetName();  pname.Append("_px");
01111    }
01112    TH3D *h1 = new TH3D(pname,GetTitle(),nx,fXaxis.GetXmin(),fXaxis.GetXmax(),ny,fYaxis.GetXmin(),fYaxis.GetXmax(),nz,fZaxis.GetXmin(),fZaxis.GetXmax());
01113    Bool_t computeErrors = kFALSE;
01114    Bool_t cequalErrors  = kFALSE;
01115    Bool_t binEntries    = kFALSE;
01116    if (opt.Contains("b")) binEntries = kTRUE;
01117    if (opt.Contains("e")) computeErrors = kTRUE;
01118    if (opt.Contains("c=e")) {cequalErrors = kTRUE; computeErrors=kFALSE;}
01119    if (computeErrors) h1->Sumw2();
01120 
01121    // Fill the projected histogram
01122    Int_t bin,binx,biny,binz;
01123    Double_t cont,err;
01124    for (binx =0;binx<=nx+1;binx++) {
01125       for (biny =0;biny<=ny+1;biny++) {
01126          for (binz =0;binz<=nz+1;binz++) {
01127             bin = GetBin(binx,biny,binz);
01128             if (binEntries)    cont = GetBinEntries(bin);
01129             else               cont = GetBinContent(bin);
01130             err   = GetBinError(bin);
01131             if (cequalErrors)  h1->SetBinContent(binx,biny,binz, err);
01132             else               h1->SetBinContent(binx,biny,binz,cont);
01133             if (computeErrors) h1->SetBinError(binx,biny,binz,err);
01134          }
01135       }
01136    }
01137    h1->SetEntries(fEntries);
01138    return h1;
01139 }
01140 
01141 //______________________________________________________________________________
01142 void TProfile3D::PutStats(Double_t *stats)
01143 {
01144    // Replace current statistics with the values in array stats
01145 
01146    TH3::PutStats(stats);
01147    fTsumwt  = stats[11];
01148    fTsumwt2 = stats[12];
01149 }
01150 
01151 //______________________________________________________________________________
01152 void TProfile3D::Reset(Option_t *option)
01153 {
01154 //*-*-*-*-*-*-*-*-*-*Reset contents of a Profile3D histogram*-*-*-*-*-*-*-*
01155 //*-*                =======================================
01156    TH3D::Reset(option);
01157    fBinSumw2.Reset();
01158    fBinEntries.Reset();
01159    TString opt = option;
01160    opt.ToUpper();
01161    if (opt.Contains("ICE")) return;
01162    fTsumwt = fTsumwt2 = 0;
01163 }
01164 
01165 //______________________________________________________________________________
01166 void TProfile3D::RebinAxis(Double_t x, TAxis *axis)
01167 {
01168 // Profile histogram is resized along axis such that x is in the axis range.
01169 // The new axis limits are recomputed by doubling iteratively
01170 // the current axis range until the specified value x is within the limits.
01171 // The algorithm makes a copy of the histogram, then loops on all bins
01172 // of the old histogram to fill the rebinned histogram.
01173 // Takes into account errors (Sumw2) if any.
01174 // The bit kCanRebin must be set before invoking this function.
01175 //  Ex:  h->SetBit(TH1::kCanRebin);
01176 
01177    TProfile3D* hold = TProfileHelper::RebinAxis(this, x, axis);
01178    if ( hold ) {
01179       fTsumwt  = hold->fTsumwt;
01180       fTsumwt2 = hold->fTsumwt2;
01181       delete hold;
01182    }
01183 }
01184 
01185 //______________________________________________________________________________
01186 void TProfile3D::SavePrimitive(ostream &out, Option_t *option /*= ""*/)
01187 {
01188    // Save primitive as a C++ statement(s) on output stream out
01189 
01190    //Note the following restrictions in the code generated:
01191    // - variable bin size not implemented
01192    // - SetErrorOption not implemented
01193 
01194 
01195    char quote = '"';
01196    out <<"   "<<endl;
01197    out <<"   "<<ClassName()<<" *";
01198 
01199    out << GetName() << " = new " << ClassName() << "(" << quote
01200        << GetName() << quote << "," << quote<< GetTitle() << quote
01201        << "," << GetXaxis()->GetNbins();
01202    out << "," << GetXaxis()->GetXmin()
01203        << "," << GetXaxis()->GetXmax();
01204    out << "," << GetYaxis()->GetNbins();
01205    out << "," << GetYaxis()->GetXmin()
01206        << "," << GetYaxis()->GetXmax();
01207    out << "," << GetZaxis()->GetNbins();
01208    out << "," << GetZaxis()->GetXmin()
01209        << "," << GetZaxis()->GetXmax();
01210    out << "," << fTmin
01211        << "," << fTmax;
01212    out << ");" << endl;
01213 
01214 
01215    // save bin entries
01216    Int_t bin;
01217    for (bin=0;bin<fNcells;bin++) {
01218       Double_t bi = GetBinEntries(bin);
01219       if (bi) {
01220          out<<"   "<<GetName()<<"->SetBinEntries("<<bin<<","<<bi<<");"<<endl;
01221       }
01222    }
01223    //save bin contents
01224    for (bin=0;bin<fNcells;bin++) {
01225       Double_t bc = fArray[bin];
01226       if (bc) {
01227          out<<"   "<<GetName()<<"->SetBinContent("<<bin<<","<<bc<<");"<<endl;
01228       }
01229    }
01230    // save bin errors
01231    if (fSumw2.fN) {
01232       for (bin=0;bin<fNcells;bin++) {
01233          Double_t be = TMath::Sqrt(fSumw2.fArray[bin]);
01234          if (be) {
01235             out<<"   "<<GetName()<<"->SetBinError("<<bin<<","<<be<<");"<<endl;
01236          }
01237       }
01238    }
01239 
01240    TH1::SavePrimitiveHelp(out, GetName(), option);
01241 }
01242 
01243 //______________________________________________________________________________
01244 void TProfile3D::Scale(Double_t c1, Option_t *option)
01245 {
01246 // *-*-*-*-*Multiply this profile2D by a constant c1*-*-*-*-*-*-*-*-*
01247 // *-*      ========================================
01248 //
01249 //   this = c1*this
01250 //
01251 // This function uses the services of TProfile3D::Add
01252 //
01253 
01254    TProfileHelper::Scale(this, c1, option);
01255 }
01256 
01257 //______________________________________________________________________________
01258 void TProfile3D::SetBinEntries(Int_t bin, Double_t w)
01259 {
01260 //*-*-*-*-*-*-*-*-*Set the number of entries in bin*-*-*-*-*-*-*-*-*-*-*-*
01261 //*-*              ================================
01262 
01263    if (bin < 0 || bin >= fNcells) return;
01264    fBinEntries.fArray[bin] = w;
01265 }
01266 
01267 //______________________________________________________________________________
01268 void TProfile3D::SetBins(Int_t nx, Double_t xmin, Double_t xmax, Int_t ny, Double_t ymin, Double_t ymax, Int_t nz, Double_t zmin, Double_t zmax)
01269 {
01270 //*-*-*-*-*-*-*-*-*Redefine  x axis parameters*-*-*-*-*-*-*-*-*-*-*-*
01271 //*-*              ===========================
01272 
01273    fXaxis.Set(nx,xmin,xmax);
01274    fYaxis.Set(ny,ymin,ymax);
01275    fZaxis.Set(ny,zmin,zmax);
01276    fNcells = (nx+2)*(ny+2)*(nz+2);
01277    fBinEntries.Set(fNcells);
01278    fSumw2.Set(fNcells);
01279    if (fBinSumw2.fN) fBinSumw2.Set(fNcells);
01280 }
01281 
01282 
01283 //______________________________________________________________________________
01284 void TProfile3D::SetBuffer(Int_t buffersize, Option_t *)
01285 {
01286 // set the buffer size in units of 8 bytes (double)
01287 
01288    if (fBuffer) {
01289       BufferEmpty();
01290       delete [] fBuffer;
01291       fBuffer = 0;
01292    }
01293    if (buffersize <= 0) {
01294       fBufferSize = 0;
01295       return;
01296    }
01297    if (buffersize < 100) buffersize = 100;
01298    fBufferSize = 1 + 5*buffersize; 
01299    fBuffer = new Double_t[fBufferSize];
01300    memset(fBuffer,0,8*fBufferSize);
01301 }
01302 
01303 //______________________________________________________________________________
01304 void TProfile3D::SetErrorOption(Option_t *option)
01305 {
01306 //*-*-*-*-*-*-*-*-*-*Set option to compute profile2D errors*-*-*-*-*-*-*-*
01307 //*-*                =======================================
01308 //
01309 //    The computation of errors is based on the parameter option:
01310 //    option:
01311 //     ' '  (Default) Errors are Spread/SQRT(N) for Spread.ne.0. ,
01312 //                      "     "  SQRT(T)/SQRT(N) for Spread.eq.0,N.gt.0 ,
01313 //                      "     "  0.  for N.eq.0
01314 //     's'            Errors are Spread  for Spread.ne.0. ,
01315 //                      "     "  SQRT(T)  for Spread.eq.0,N.gt.0 ,
01316 //                      "     "  0.  for N.eq.0
01317 //     'i'            Errors are Spread/SQRT(N) for Spread.ne.0. ,
01318 //                      "     "  1./SQRT(12.*N) for Spread.eq.0,N.gt.0 ,
01319 //                      "     "  0.  for N.eq.0
01320 //   See TProfile3D::BuildOptions for explanation of all options
01321 
01322    TString opt = option;
01323    opt.ToLower();
01324    fErrorMode = kERRORMEAN;
01325    if (opt.Contains("s")) fErrorMode = kERRORSPREAD;
01326    if (opt.Contains("i")) fErrorMode = kERRORSPREADI;
01327    if (opt.Contains("g")) fErrorMode = kERRORSPREADG;
01328 }
01329 
01330 //______________________________________________________________________________
01331 void TProfile3D::Sumw2()
01332 {
01333    // Create structure to store sum of squares of weights per bin  *-*-*-*-*-*-*-*
01334    //   This is needed to compute  the correct statistical quantities  
01335    //    of a profile filled with weights 
01336    //  
01337    //
01338    //  This function is automatically called when the histogram is created
01339    //  if the static function TH1::SetDefaultSumw2 has been called before.
01340 
01341    TProfileHelper::Sumw2(this);
01342 }

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