Logo ROOT   6.10/09
Reference Guide
TProfile2D.cxx
Go to the documentation of this file.
1 // @(#)root/hist:$Id$
2 // Author: Rene Brun 16/04/2000
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 #include "TProfile2D.h"
13 #include "TBuffer.h"
14 #include "TMath.h"
15 #include "THLimitsFinder.h"
16 #include "Riostream.h"
17 #include "TVirtualPad.h"
18 #include "TError.h"
19 #include "TClass.h"
20 
21 #include "TProfileHelper.h"
22 
24 
26 
27 /** \class TProfile2D
28  \ingroup Hist
29  Profile2D histograms are used to display the mean
30  value of Z and its RMS for each cell in X,Y.
31  Profile2D histograms are in many cases an
32  elegant replacement of three-dimensional histograms : the inter-relation of three
33  measured quantities X, Y and Z can always be visualized by a three-dimensional
34  histogram or scatter-plot; its representation on the line-printer is not particularly
35  satisfactory, except for sparse data. If Z is an unknown (but single-valued)
36  approximate function of X,Y this function is displayed by a profile2D histogram with
37  much better precision than by a scatter-plot.
38 
39  The following formulae show the cumulated contents (capital letters) and the values
40  displayed by the printing or plotting routines (small letters) of the elements for cell I, J.
41 
42  2
43  H(I,J) = sum Z E(I,J) = sum Z
44  l(I,J) = sum l L(I,J) = sum l
45  h(I,J) = H(I,J)/L(I,J) s(I,J) = sqrt(E(I,J)/L(I,J)- h(I,J)**2)
46  e(I,J) = s(I,J)/sqrt(L(I,J))
47 
48  In the special case where s(I,J) is zero (eg, case of 1 entry only in one cell)
49  the bin error e(I,J) is computed from the average of the s(I,J) for all cells
50  if the static function TProfile2D::Approximate has been called.
51  This simple/crude approximation was suggested in order to keep the cell
52  during a fit operation. But note that this approximation is not the default behaviour.
53 
54  Example of a profile2D histogram
55  ~~~~{.cpp}
56  {
57  auto c1 = new TCanvas("c1","Profile histogram example",200,10,700,500);
58  auto hprof2d = new TProfile2D("hprof2d","Profile of pz versus px and py",40,-4,4,40,-4,4,0,20);
59  Float_t px, py, pz;
60  for ( Int_t i=0; i<25000; i++) {
61  gRandom->Rannor(px,py);
62  pz = px*px + py*py;
63  hprof2d->Fill(px,py,pz,1);
64  }
65  hprof2d->Draw();
66  }
67  ~~~~
68 */
69 
70 ////////////////////////////////////////////////////////////////////////////////
71 /// Default constructor for Profile2D histograms.
72 
74 {
75  fTsumwz = fTsumwz2 = 0;
76  fScaling = kFALSE;
77  BuildOptions(0,0,"");
78 }
79 
80 ////////////////////////////////////////////////////////////////////////////////
81 /// Default destructor for Profile2D histograms.
82 
84 {
85 }
86 
87 ////////////////////////////////////////////////////////////////////////////////
88 /// Normal Constructor for Profile histograms.
89 ///
90 /// The first eight parameters are similar to TH2D::TH2D.
91 /// All values of z are accepted at filling time.
92 /// To fill a profile2D histogram, one must use TProfile2D::Fill function.
93 ///
94 /// Note that when filling the profile histogram the function Fill
95 /// checks if the variable z is between fZmin and fZmax.
96 /// If a minimum or maximum value is set for the Z scale before filling,
97 /// then all values below zmin or above zmax will be discarded.
98 /// Setting the minimum or maximum value for the Z scale before filling
99 /// has the same effect as calling the special TProfile2D constructor below
100 /// where zmin and zmax are specified.
101 ///
102 /// H(I,J) is printed as the cell contents. The errors computed are s(I,J) if CHOPT='S'
103 /// (spread option), or e(I,J) if CHOPT=' ' (error on mean).
104 ///
105 /// See TProfile2D::BuildOptions for explanation of parameters
106 ///
107 /// see other constructors below with all possible combinations of
108 /// fix and variable bin size like in TH2D.
109 
110 TProfile2D::TProfile2D(const char *name,const char *title,Int_t nx,Double_t xlow,Double_t xup,Int_t ny,Double_t ylow,Double_t yup,Option_t *option)
111 : TH2D(name,title,nx,xlow,xup,ny,ylow,yup)
112 {
113  BuildOptions(0,0,option);
114  if (xlow >= xup || ylow >= yup) SetBuffer(fgBufferSize);
115 }
116 
117 ////////////////////////////////////////////////////////////////////////////////
118 /// Create a 2-D Profile with variable bins in X and fix bins in Y.
119 
120 TProfile2D::TProfile2D(const char *name,const char *title,Int_t nx,const Double_t *xbins,Int_t ny,Double_t ylow,Double_t yup,Option_t *option)
121 : TH2D(name,title,nx,xbins,ny,ylow,yup)
122 {
123  BuildOptions(0,0,option);
124 }
125 
126 ////////////////////////////////////////////////////////////////////////////////
127 /// Create a 2-D Profile with fix bins in X and variable bins in Y.
128 
129 TProfile2D::TProfile2D(const char *name,const char *title,Int_t nx,Double_t xlow,Double_t xup,Int_t ny,const Double_t *ybins,Option_t *option)
130 : TH2D(name,title,nx,xlow,xup,ny,ybins)
131 {
132  BuildOptions(0,0,option);
133 }
134 
135 ////////////////////////////////////////////////////////////////////////////////
136 /// Create a 2-D Profile with variable bins in X and variable bins in Y.
137 
138 TProfile2D::TProfile2D(const char *name,const char *title,Int_t nx,const Double_t *xbins,Int_t ny,const Double_t *ybins,Option_t *option)
139 : TH2D(name,title,nx,xbins,ny,ybins)
140 {
141  BuildOptions(0,0,option);
142 }
143 
144 ////////////////////////////////////////////////////////////////////////////////
145 /// Constructor for Profile2D histograms with range in z.
146 ///
147 /// The first eight parameters are similar to TH2D::TH2D.
148 /// Only the values of Z between ZMIN and ZMAX will be considered at filling time.
149 /// zmin and zmax will also be the maximum and minimum values
150 /// on the z scale when drawing the profile2D.
151 ///
152 /// See TProfile2D::BuildOptions for more explanations on errors
153 
154 TProfile2D::TProfile2D(const char *name,const char *title,Int_t nx,Double_t xlow,Double_t xup,Int_t ny, Double_t ylow,Double_t yup,Double_t zlow,Double_t zup,Option_t *option)
155 : TH2D(name,title,nx,xlow,xup,ny,ylow,yup)
156 {
157  BuildOptions(zlow,zup,option);
158  if (xlow >= xup || ylow >= yup) SetBuffer(fgBufferSize);
159 }
160 
161 
162 ////////////////////////////////////////////////////////////////////////////////
163 /// Set Profile2D histogram structure and options.
164 ///
165 /// - zmin: minimum value allowed for z
166 /// - zmax: maximum value allowed for z
167 /// if (zmin = zmax = 0) there are no limits on the allowed z values (zmin = -inf, zmax = +inf)
168 ///
169 /// - option: this is the option for the computation of the t error of the profile ( TProfile2D::GetBinError )
170 /// possible values for the options are documented in TProfile2D::SetErrorOption
171 ///
172 /// See TProfile::BuildOptions for a detailed description
173 
175 {
176 
177  SetErrorOption(option);
178 
179  // create extra profile data structure (bin entries/ y^2 and sum of weight square)
181 
182  fZmin = zmin;
183  fZmax = zmax;
184  fScaling = kFALSE;
185  fTsumwz = fTsumwz2 = 0;
186 }
187 
188 ////////////////////////////////////////////////////////////////////////////////
189 /// Copy constructor.
190 
192 {
193  ((TProfile2D&)profile).Copy(*this);
194 }
195 
196 ////////////////////////////////////////////////////////////////////////////////
197 /// Performs the operation: `this = this + c1*f1` .
198 
200 {
201  Error("Add","Function not implemented for TProfile2D");
202  return kFALSE;
203 }
204 
205 ////////////////////////////////////////////////////////////////////////////////
206 /// Performs the operation: `this = this + c1*h1` .
207 
209 {
210  if (!h1) {
211  Error("Add","Attempt to add a non-existing profile");
212  return kFALSE;
213  }
214  if (!h1->InheritsFrom(TProfile2D::Class())) {
215  Error("Add","Attempt to add a non-profile2D object");
216  return kFALSE;
217  }
218 
219  return TProfileHelper::Add(this, this, h1, 1, c1);
220 }
221 
222 ////////////////////////////////////////////////////////////////////////////////
223 /// Replace contents of this profile2D by the addition of h1 and h2.
224 ///
225 /// `this = c1*h1 + c2*h2`
226 
228 {
229  if (!h1 || !h2) {
230  Error("Add","Attempt to add a non-existing profile");
231  return kFALSE;
232  }
233  if (!h1->InheritsFrom(TProfile2D::Class())) {
234  Error("Add","Attempt to add a non-profile2D object");
235  return kFALSE;
236  }
237  if (!h2->InheritsFrom(TProfile2D::Class())) {
238  Error("Add","Attempt to add a non-profile2D object");
239  return kFALSE;
240  }
241  return TProfileHelper::Add(this, h1, h2, c1, c2);
242 }
243 
244 ////////////////////////////////////////////////////////////////////////////////
245 /// Static function, set the fgApproximate flag.
246 ///
247 /// When the flag is true, the function GetBinError
248 /// will approximate the bin error with the average profile error on all bins
249 /// in the following situation only
250 /// - the number of bins in the profile2D is less than 10404 (eg 100x100)
251 /// - the bin number of entries is small ( <5)
252 /// - the estimated bin error is extremely small compared to the bin content
253 /// (see TProfile2D::GetBinError)
254 
256 {
257  fgApproximate = approx;
258 }
259 
260 ////////////////////////////////////////////////////////////////////////////////
261 /// Fill histogram with all entries in the buffer.
262 ///
263 /// - action = -1 histogram is reset and refilled from the buffer (called by THistPainter::Paint)
264 /// - action = 0 histogram is filled from the buffer
265 /// - action = 1 histogram is filled and buffer is deleted
266 /// The buffer is automatically deleted when the number of entries
267 /// in the buffer is greater than the number of entries in the histogram
268 
270 {
271  // do we need to compute the bin size?
272  if (!fBuffer) return 0;
273  Int_t nbentries = (Int_t)fBuffer[0];
274  if (!nbentries) return 0;
275  Double_t *buffer = fBuffer;
276  if (nbentries < 0) {
277  if (action == 0) return 0;
278  nbentries = -nbentries;
279  fBuffer=0;
280  Reset("ICES"); // reset without deleting the functions
281  fBuffer = buffer;
282  }
284  //find min, max of entries in buffer
285  Double_t xmin = fBuffer[2];
286  Double_t xmax = xmin;
287  Double_t ymin = fBuffer[3];
288  Double_t ymax = ymin;
289  for (Int_t i=1;i<nbentries;i++) {
290  Double_t x = fBuffer[4*i+2];
291  if (x < xmin) xmin = x;
292  if (x > xmax) xmax = x;
293  Double_t y = fBuffer[4*i+3];
294  if (y < ymin) ymin = y;
295  if (y > ymax) ymax = y;
296  }
297  if (fXaxis.GetXmax() <= fXaxis.GetXmin() || fYaxis.GetXmax() <= fYaxis.GetXmin()) {
298  THLimitsFinder::GetLimitsFinder()->FindGoodLimits(this,xmin,xmax,ymin,ymax);
299  } else {
300  fBuffer = 0;
301  Int_t keep = fBufferSize; fBufferSize = 0;
302  if (xmin < fXaxis.GetXmin()) ExtendAxis(xmin,&fXaxis);
303  if (xmax >= fXaxis.GetXmax()) ExtendAxis(xmax,&fXaxis);
304  if (ymin < fYaxis.GetXmin()) ExtendAxis(ymin,&fYaxis);
305  if (ymax >= fYaxis.GetXmax()) ExtendAxis(ymax,&fYaxis);
306  fBuffer = buffer;
307  fBufferSize = keep;
308  }
309  }
310 
311  fBuffer = 0;
312  for (Int_t i=0;i<nbentries;i++) {
313  Fill(buffer[4*i+2],buffer[4*i+3],buffer[4*i+4],buffer[4*i+1]);
314  }
315  fBuffer = buffer;
316 
317  if (action > 0) { delete [] fBuffer; fBuffer = 0; fBufferSize = 0;}
318  else {
319  if (nbentries == (Int_t)fEntries) fBuffer[0] = -nbentries;
320  else fBuffer[0] = 0;
321  }
322  return nbentries;
323 }
324 
325 ////////////////////////////////////////////////////////////////////////////////
326 /// Accumulate arguments in buffer.
327 ///
328 /// When buffer is full, empty the buffer.
329 ///
330 /// - fBuffer[0] = number of entries in buffer
331 /// - fBuffer[1] = w of first entry
332 /// - fBuffer[2] = x of first entry
333 /// - fBuffer[3] = y of first entry
334 /// - fBuffer[4] = z of first entry
335 
337 {
338  if (!fBuffer) return -3;
339  Int_t nbentries = (Int_t)fBuffer[0];
340  if (nbentries < 0) {
341  nbentries = -nbentries;
342  fBuffer[0] = nbentries;
343  if (fEntries > 0) {
344  Double_t *buffer = fBuffer; fBuffer=0;
345  Reset("ICES"); // reset without deleting the functions
346  fBuffer = buffer;
347  }
348  }
349  if (4*nbentries+4 >= fBufferSize) {
350  BufferEmpty(1);
351  return Fill(x,y,z,w);
352  }
353  fBuffer[4*nbentries+1] = w;
354  fBuffer[4*nbentries+2] = x;
355  fBuffer[4*nbentries+3] = y;
356  fBuffer[4*nbentries+4] = z;
357  fBuffer[0] += 1;
358  return -2;
359 }
360 
361 ////////////////////////////////////////////////////////////////////////////////
362 /// Copy a Profile2D histogram to a new profile2D histogram.
363 
364 void TProfile2D::Copy(TObject &obj) const
365 {
366  try {
367  TProfile2D & pobj = dynamic_cast<TProfile2D&>(obj);
368 
369  TH2D::Copy(pobj);
371  fBinSumw2.Copy(pobj.fBinSumw2);
372  for (int bin=0;bin<fNcells;bin++) {
373  pobj.fArray[bin] = fArray[bin];
374  pobj.fSumw2.fArray[bin] = fSumw2.fArray[bin];
375  }
376  pobj.fZmin = fZmin;
377  pobj.fZmax = fZmax;
378  pobj.fScaling = fScaling;
379  pobj.fErrorMode = fErrorMode;
380  pobj.fTsumwz = fTsumwz;
381  pobj.fTsumwz2 = fTsumwz2;
382 
383  } catch(...) {
384  Fatal("Copy","Cannot copy a TProfile2D in a %s",obj.IsA()->GetName());
385  }
386 
387 }
388 
389 ////////////////////////////////////////////////////////////////////////////////
390 /// Performs the operation: `this = this/(c1*f1)` .
391 /// This function is not implemented
392 
394 {
395  Error("Divide","Function not implemented for TProfile2D");
396  return kFALSE;
397 }
398 
399 ////////////////////////////////////////////////////////////////////////////////
400 /// Divide this profile2D by h1.
401 ///
402 /// `this = this/h1`
403 ///
404 ///This function return kFALSE if the divide operation failed
405 
407 {
408 
409  if (!h1) {
410  Error("Divide","Attempt to divide a non-existing profile2D");
411  return kFALSE;
412  }
413  if (!h1->InheritsFrom(TProfile2D::Class())) {
414  Error("Divide","Attempt to divide a non-profile2D object");
415  return kFALSE;
416  }
417  TProfile2D *p1 = (TProfile2D*)h1;
418 
419  // delete buffer if it is there since it will become invalid
420  if (fBuffer) BufferEmpty(1);
421 
422  // Check profile compatibility
423  Int_t nx = GetNbinsX();
424  if (nx != p1->GetNbinsX()) {
425  Error("Divide","Attempt to divide profiles with different number of bins");
426  return kFALSE;
427  }
428  Int_t ny = GetNbinsY();
429  if (ny != p1->GetNbinsY()) {
430  Error("Divide","Attempt to divide profiles with different number of bins");
431  return kFALSE;
432  }
433 
434  // Reset statistics
435  fEntries = fTsumw = fTsumw2 = fTsumwx = fTsumwx2 = 0;
436 
437  // Loop on bins (including underflows/overflows)
438  Int_t bin,binx,biny;
439  Double_t *cu1 = p1->GetW();
440  Double_t *er1 = p1->GetW2();
441  Double_t *en1 = p1->GetB();
442  Double_t c0,c1,w,z,x,y;
443  for (binx =0;binx<=nx+1;binx++) {
444  for (biny =0;biny<=ny+1;biny++) {
445  bin = biny*(fXaxis.GetNbins()+2) + binx;
446  c0 = fArray[bin];
447  c1 = cu1[bin];
448  if (c1) w = c0/c1;
449  else w = 0;
450  fArray[bin] = w;
451  z = TMath::Abs(w);
452  x = fXaxis.GetBinCenter(binx);
453  y = fYaxis.GetBinCenter(biny);
454  fEntries++;
455  fTsumw += z;
456  fTsumw2 += z*z;
457  fTsumwx += z*x;
458  fTsumwx2 += z*x*x;
459  fTsumwy += z*y;
460  fTsumwy2 += z*y*y;
461  fTsumwxy += z*x*y;
462  fTsumwz += z;
463  fTsumwz2 += z*z;
464  Double_t e0 = fSumw2.fArray[bin];
465  Double_t e1 = er1[bin];
466  Double_t c12= c1*c1;
467  if (!c1) fSumw2.fArray[bin] = 0;
468  else fSumw2.fArray[bin] = (e0*c1*c1 + e1*c0*c0)/(c12*c12);
469  if (!en1[bin]) fBinEntries.fArray[bin] = 0;
470  else fBinEntries.fArray[bin] /= en1[bin];
471  }
472  }
473  // maintaining the correct sum of weights square is not supported when dividing
474  // bin error resulting from division of profile needs to be checked
475  if (fBinSumw2.fN) {
476  Warning("Divide","Cannot preserve during the division of profiles the sum of bin weight square");
477  fBinSumw2 = TArrayD();
478  }
479  return kTRUE;
480 }
481 
482 ////////////////////////////////////////////////////////////////////////////////
483 /// Replace contents of this profile2D by the division of h1 by h2.
484 ///
485 /// `this = c1*h1/(c2*h2)`
486 ///
487 /// This function return kFALSE if the divide operation failed
488 
490 {
491  TString opt = option;
492  opt.ToLower();
493  Bool_t binomial = kFALSE;
494  if (opt.Contains("b")) binomial = kTRUE;
495  if (!h1 || !h2) {
496  Error("Divide","Attempt to divide a non-existing profile2D");
497  return kFALSE;
498  }
499  if (!h1->InheritsFrom(TProfile2D::Class())) {
500  Error("Divide","Attempt to divide a non-profile2D object");
501  return kFALSE;
502  }
503  TProfile2D *p1 = (TProfile2D*)h1;
504  if (!h2->InheritsFrom(TProfile2D::Class())) {
505  Error("Divide","Attempt to divide a non-profile2D object");
506  return kFALSE;
507  }
508  TProfile2D *p2 = (TProfile2D*)h2;
509 
510  // delete buffer if it is there since it will become invalid
511  if (fBuffer) BufferEmpty(1);
512 
513  // Check histogram compatibility
514  Int_t nx = GetNbinsX();
515  if (nx != p1->GetNbinsX() || nx != p2->GetNbinsX()) {
516  Error("Divide","Attempt to divide profiles with different number of bins");
517  return kFALSE;
518  }
519  Int_t ny = GetNbinsY();
520  if (ny != p1->GetNbinsY() || ny != p2->GetNbinsY()) {
521  Error("Divide","Attempt to divide profiles with different number of bins");
522  return kFALSE;
523  }
524  if (!c2) {
525  Error("Divide","Coefficient of dividing profile cannot be zero");
526  return kFALSE;
527  }
528 
529  // Reset statistics
530  fEntries = fTsumw = fTsumw2 = fTsumwx = fTsumwx2 = 0;
531 
532  // Loop on bins (including underflows/overflows)
533  Int_t bin,binx,biny;
534  Double_t *cu1 = p1->GetW();
535  Double_t *cu2 = p2->GetW();
536  Double_t *er1 = p1->GetW2();
537  Double_t *er2 = p2->GetW2();
538  Double_t *en1 = p1->GetB();
539  Double_t *en2 = p2->GetB();
540  Double_t b1,b2,w,z,x,y,ac1,ac2;
541  ac1 = TMath::Abs(c1);
542  ac2 = TMath::Abs(c2);
543  for (binx =0;binx<=nx+1;binx++) {
544  for (biny =0;biny<=ny+1;biny++) {
545  bin = biny*(fXaxis.GetNbins()+2) + binx;
546  b1 = cu1[bin];
547  b2 = cu2[bin];
548  if (b2) w = c1*b1/(c2*b2);
549  else w = 0;
550  fArray[bin] = w;
551  z = TMath::Abs(w);
552  x = fXaxis.GetBinCenter(binx);
553  y = fYaxis.GetBinCenter(biny);
554  fEntries++;
555  fTsumw += z;
556  fTsumw2 += z*z;
557  fTsumwx += z*x;
558  fTsumwx2 += z*x*x;
559  fTsumwy += z*y;
560  fTsumwy2 += z*y*y;
561  fTsumwxy += z*x*y;
562  fTsumwz += z;
563  fTsumwz2 += z*z;
564  Double_t e1 = er1[bin];
565  Double_t e2 = er2[bin];
566  //Double_t b22= b2*b2*d2;
567  Double_t b22= b2*b2*TMath::Abs(c2);
568  if (!b2) fSumw2.fArray[bin] = 0;
569  else {
570  if (binomial) {
571  fSumw2.fArray[bin] = TMath::Abs(w*(1-w)/(c2*b2));
572  } else {
573  fSumw2.fArray[bin] = ac1*ac2*(e1*b2*b2 + e2*b1*b1)/(b22*b22);
574  }
575  }
576  if (!en2[bin]) fBinEntries.fArray[bin] = 0;
577  else fBinEntries.fArray[bin] = en1[bin]/en2[bin];
578  }
579  }
580  return kTRUE;
581 }
582 
583 ////////////////////////////////////////////////////////////////////////////////
584 /// Fill a Profile2D histogram (no weights).
585 
587 {
588  if (fBuffer) return BufferFill(x,y,z,1);
589 
590  Int_t bin,binx,biny;
591 
592  if (fZmin != fZmax) {
593  if (z <fZmin || z> fZmax || TMath::IsNaN(z) ) return -1;
594  }
595 
596  fEntries++;
597  binx =fXaxis.FindBin(x);
598  biny =fYaxis.FindBin(y);
599  if (binx <0 || biny <0) return -1;
600  bin = GetBin(binx, biny);
601  fArray[bin] += z;
602  fSumw2.fArray[bin] += z*z;
603  fBinEntries.fArray[bin] += 1;
604  if (fBinSumw2.fN) fBinSumw2.fArray[bin] += 1;
605  if (binx == 0 || binx > fXaxis.GetNbins()) {
606  if (!fgStatOverflows) return -1;
607  }
608  if (biny == 0 || biny > fYaxis.GetNbins()) {
609  if (!fgStatOverflows) return -1;
610  }
611  ++fTsumw;
612  ++fTsumw2;
613  fTsumwx += x;
614  fTsumwx2 += x*x;
615  fTsumwy += y;
616  fTsumwy2 += y*y;
617  fTsumwxy += x*y;
618  fTsumwz += z;
619  fTsumwz2 += z*z;
620  return bin;
621 }
622 
623 ////////////////////////////////////////////////////////////////////////////////
624 /// Fill a Profile2D histogram (no weights).
625 
627 {
628  Int_t bin,binx,biny;
629 
630  if (fZmin != fZmax) {
631  if (z <fZmin || z> fZmax || TMath::IsNaN(z)) return -1;
632  }
633 
634  fEntries++;
635  binx =fXaxis.FindBin(x);
636  biny =fYaxis.FindBin(namey);
637  if (binx <0 || biny <0) return -1;
638  bin = biny*(fXaxis.GetNbins()+2) + binx;
639  AddBinContent(bin, z);
640  fSumw2.fArray[bin] += (Double_t)z*z;
641  fBinEntries.fArray[bin] += 1;
642  if (fBinSumw2.fN) fBinSumw2.fArray[bin] += 1;
643  if (binx == 0 || binx > fXaxis.GetNbins()) {
644  if (!fgStatOverflows) return -1;
645  }
646  if (biny == 0 || biny > fYaxis.GetNbins()) return -1;
647  Double_t y = fYaxis.GetBinCenter(biny);
648  ++fTsumw;
649  ++fTsumw2;
650  fTsumwx += x;
651  fTsumwx2 += x*x;
652  fTsumwy += y;
653  fTsumwy2 += y*y;
654  fTsumwxy += x*y;
655  fTsumwz += z;
656  fTsumwz2 += z*z;
657  return bin;
658 }
659 
660 ////////////////////////////////////////////////////////////////////////////////
661 /// Fill a Profile2D histogram (no weights).
662 
663 Int_t TProfile2D::Fill(const char *namex, const char *namey, Double_t z)
664 {
665  Int_t bin,binx,biny;
666 
667  if (fZmin != fZmax) {
668  if (z <fZmin || z> fZmax || TMath::IsNaN(z) ) return -1;
669  }
670 
671  fEntries++;
672  binx =fXaxis.FindBin(namex);
673  biny =fYaxis.FindBin(namey);
674  if (binx <0 || biny <0) return -1;
675  bin = biny*(fXaxis.GetNbins()+2) + binx;
676  AddBinContent(bin, z);
677  fSumw2.fArray[bin] += (Double_t)z*z;
678  fBinEntries.fArray[bin] += 1;
679  if (fBinSumw2.fN) fBinSumw2.fArray[bin] += 1;
680  if (binx == 0 || binx > fXaxis.GetNbins()) return -1;
681  if (biny == 0 || biny > fYaxis.GetNbins()) return -1;
682  Double_t x = fYaxis.GetBinCenter(binx);
683  Double_t y = fYaxis.GetBinCenter(biny);
684  ++fTsumw;
685  ++fTsumw2;
686  fTsumwx += x;
687  fTsumwx2 += x*x;
688  fTsumwy += y;
689  fTsumwy2 += y*y;
690  fTsumwxy += x*y;
691  fTsumwz += z;
692  fTsumwz2 += z*z;
693  return bin;
694 }
695 
696 ////////////////////////////////////////////////////////////////////////////////
697 /// Fill a Profile2D histogram (no weights).
698 
700 {
701  Int_t bin,binx,biny;
702 
703  if (fZmin != fZmax) {
704  if (z <fZmin || z> fZmax || TMath::IsNaN(z)) return -1;
705  }
706 
707  fEntries++;
708  binx =fXaxis.FindBin(namex);
709  biny =fYaxis.FindBin(y);
710  if (binx <0 || biny <0) return -1;
711  bin = biny*(fXaxis.GetNbins()+2) + binx;
712  AddBinContent(bin, z);
713  fSumw2.fArray[bin] += (Double_t)z*z;
714  fBinEntries.fArray[bin] += 1;
715  if (fBinSumw2.fN) fBinSumw2.fArray[bin] += 1;
716  if (binx == 0 || binx > fXaxis.GetNbins()) return -1;
717  if (biny == 0 || biny > fYaxis.GetNbins()) {
718  if (!fgStatOverflows) return -1;
719  }
720  Double_t x = fYaxis.GetBinCenter(binx);
721  ++fTsumw;
722  ++fTsumw2;
723  fTsumwx += x;
724  fTsumwx2 += x*x;
725  fTsumwy += y;
726  fTsumwy2 += y*y;
727  fTsumwxy += x*y;
728  fTsumwz += z;
729  fTsumwz2 += z*z;
730  return bin;
731 }
732 
733 ////////////////////////////////////////////////////////////////////////////////
734 /// Fill a Profile2D histogram with weights.
735 
737 {
738  if (fBuffer) return BufferFill(x,y,z,w);
739 
740  Int_t bin,binx,biny;
741 
742  if (fZmin != fZmax) {
743  if (z <fZmin || z> fZmax || TMath::IsNaN(z)) return -1;
744  }
745 
746  Double_t u= w;
747  fEntries++;
748  binx =fXaxis.FindBin(x);
749  biny =fYaxis.FindBin(y);
750  if (binx <0 || biny <0) return -1;
751  bin = biny*(fXaxis.GetNbins()+2) + binx;
752  AddBinContent(bin, u*z);
753  fSumw2.fArray[bin] += u*z*z;
754  if (!fBinSumw2.fN && u != 1.0 && !TestBit(TH1::kIsNotW)) Sumw2(); // must be called before accumulating the entries
755  if (fBinSumw2.fN) fBinSumw2.fArray[bin] += u*u;
756  fBinEntries.fArray[bin] += u;
757  if (binx == 0 || binx > fXaxis.GetNbins()) {
758  if (!fgStatOverflows) return -1;
759  }
760  if (biny == 0 || biny > fYaxis.GetNbins()) {
761  if (!fgStatOverflows) return -1;
762  }
763  fTsumw += u;
764  fTsumw2 += u*u;
765  fTsumwx += u*x;
766  fTsumwx2 += u*x*x;
767  fTsumwy += u*y;
768  fTsumwy2 += u*y*y;
769  fTsumwxy += u*x*y;
770  fTsumwz += u*z;
771  fTsumwz2 += u*z*z;
772  return bin;
773 }
774 
775 ////////////////////////////////////////////////////////////////////////////////
776 /// Return bin content of a Profile2D histogram.
777 
779 {
780  if (fBuffer) ((TProfile2D*)this)->BufferEmpty();
781 
782  if (bin < 0 || bin >= fNcells) return 0;
783  if (fBinEntries.fArray[bin] == 0) return 0;
784  if (!fArray) return 0;
785  return fArray[bin]/fBinEntries.fArray[bin];
786 }
787 
788 ////////////////////////////////////////////////////////////////////////////////
789 /// Return bin entries of a Profile2D histogram.
790 
792 {
793  if (fBuffer) ((TProfile2D*)this)->BufferEmpty();
794 
795  if (bin < 0 || bin >= fNcells) return 0;
796  return fBinEntries.fArray[bin];
797 }
798 
799 ////////////////////////////////////////////////////////////////////////////////
800 /// Return bin effective entries for a weighted filled Profile histogram.
801 /// In case of an unweighted profile, it is equivalent to the number of entries per bin
802 /// The effective entries is defined as the square of the sum of the weights divided by the
803 /// sum of the weights square.
804 /// TProfile::Sumw2() must be called before filling the profile with weights.
805 /// Only by calling this method the sum of the square of the weights per bin is stored.
806 
808 {
809  return TProfileHelper::GetBinEffectiveEntries(this, bin);
810 }
811 
812 ////////////////////////////////////////////////////////////////////////////////
813 /// Return bin error of a Profile2D histogram.
814 ///
815 /// ### Computing errors: A moving field
816 ///
817 /// The computation of errors for a TProfile2D has evolved with the versions
818 /// of ROOT. The difficulty is in computing errors for bins with low statistics.
819 /// - prior to version 3.10, we had no special treatment of low statistic bins.
820 /// As a result, these bins had huge errors. The reason is that the
821 /// expression eprim2 is very close to 0 (rounding problems) or 0.
822 /// - The algorithm is modified/protected for the case
823 /// when a TProfile2D is projected (ProjectionX). The previous algorithm
824 /// generated a N^2 problem when projecting a TProfile2D with a large number of
825 /// bins (eg 100000).
826 /// - in version 3.10/02, a new static function TProfile::Approximate
827 /// is introduced to enable or disable (default) the approximation.
828 /// (see also comments in TProfile::GetBinError)
829 
831 {
832  return TProfileHelper::GetBinError((TProfile2D*)this, bin);
833 }
834 
835 ////////////////////////////////////////////////////////////////////////////////
836 /// Return option to compute profile2D errors.
837 
839 {
840  if (fErrorMode == kERRORSPREAD) return "s";
841  if (fErrorMode == kERRORSPREADI) return "i";
842  if (fErrorMode == kERRORSPREADG) return "g";
843  return "";
844 }
845 
846 ////////////////////////////////////////////////////////////////////////////////
847 /// Fill the array stats from the contents of this profile.
848 /// The array stats must be correctly dimensioned in the calling program.
849 ///
850 /// - stats[0] = sumw
851 /// - stats[1] = sumw2
852 /// - stats[2] = sumwx
853 /// - stats[3] = sumwx2
854 /// - stats[4] = sumwy
855 /// - stats[5] = sumwy2
856 /// - stats[6] = sumwxy
857 /// - stats[7] = sumwz
858 /// - stats[8] = sumwz2
859 ///
860 /// If no axis-subrange is specified (via TAxis::SetRange), the array stats
861 /// is simply a copy of the statistics quantities computed at filling time.
862 /// If a sub-range is specified, the function recomputes these quantities
863 /// from the bin contents in the current axis range.
864 
865 void TProfile2D::GetStats(Double_t *stats) const
866 {
867  if (fBuffer) ((TProfile2D*)this)->BufferEmpty();
868 
869  // Loop on bins
871  Int_t bin, binx, biny;
872  Double_t w, w2;
873  Double_t x,y;
874  for (bin=0;bin<9;bin++) stats[bin] = 0;
875  if (!fBinEntries.fArray) return;
876  Int_t firstBinX = fXaxis.GetFirst();
877  Int_t lastBinX = fXaxis.GetLast();
878  Int_t firstBinY = fYaxis.GetFirst();
879  Int_t lastBinY = fYaxis.GetLast();
880  // include underflow/overflow if TH1::StatOverflows(kTRUE) in case no range is set on the axis
881  if (fgStatOverflows) {
882  if ( !fXaxis.TestBit(TAxis::kAxisRange) ) {
883  if (firstBinX == 1) firstBinX = 0;
884  if (lastBinX == fXaxis.GetNbins() ) lastBinX += 1;
885  }
886  if ( !fYaxis.TestBit(TAxis::kAxisRange) ) {
887  if (firstBinY == 1) firstBinY = 0;
888  if (lastBinY == fYaxis.GetNbins() ) lastBinY += 1;
889  }
890  }
891  for (biny = firstBinY; biny <= lastBinY; biny++) {
892  y = fYaxis.GetBinCenter(biny);
893  for (binx = firstBinX; binx <= lastBinX; binx++) {
894  bin = GetBin(binx,biny);
895  w = fBinEntries.fArray[bin];
896  w2 = (fBinSumw2.fN ? fBinSumw2.fArray[bin] : w );
897  x = fXaxis.GetBinCenter(binx);
898  stats[0] += w;
899  stats[1] += w2;
900  stats[2] += w*x;
901  stats[3] += w*x*x;
902  stats[4] += w*y;
903  stats[5] += w*y*y;
904  stats[6] += w*x*y;
905  stats[7] += fArray[bin];
906  stats[8] += fSumw2.fArray[bin];
907  }
908  }
909  } else {
910  stats[0] = fTsumw;
911  stats[1] = fTsumw2;
912  stats[2] = fTsumwx;
913  stats[3] = fTsumwx2;
914  stats[4] = fTsumwy;
915  stats[5] = fTsumwy2;
916  stats[6] = fTsumwxy;
917  stats[7] = fTsumwz;
918  stats[8] = fTsumwz2;
919  }
920 }
921 
922 ////////////////////////////////////////////////////////////////////////////////
923 /// Reduce the number of bins for this axis to the number of bins having a label.
924 
926 {
928 }
929 
930 ////////////////////////////////////////////////////////////////////////////////
931 /// Double the number of bins for axis.
932 /// Refill histogram
933 /// This function is called by TAxis::FindBin(const char *label)
934 
936 {
938 }
939 
940 ////////////////////////////////////////////////////////////////////////////////
941 /// Set option(s) to draw axis with labels.
942 ///
943 /// option might have the following values:
944 ///
945 /// - "a" sort by alphabetic order
946 /// - ">" sort by decreasing values
947 /// - "<" sort by increasing values
948 /// - "h" draw labels horizontal
949 /// - "v" draw labels vertical
950 /// - "u" draw labels up (end of label right adjusted)
951 /// - "d" draw labels down (start of label left adjusted)
952 
954 {
955 
956  TAxis *axis = GetXaxis();
957  if (ax[0] == 'y' || ax[0] == 'Y') axis = GetYaxis();
958  THashList *labels = axis->GetLabels();
959  if (!labels) {
960  Warning("LabelsOption","Cannot sort. No labels");
961  return;
962  }
963  TString opt = option;
964  opt.ToLower();
965  if (opt.Contains("h")) {
966  axis->SetBit(TAxis::kLabelsHori);
969  axis->ResetBit(TAxis::kLabelsUp);
970  }
971  if (opt.Contains("v")) {
972  axis->SetBit(TAxis::kLabelsVert);
975  axis->ResetBit(TAxis::kLabelsUp);
976  }
977  if (opt.Contains("u")) {
978  axis->SetBit(TAxis::kLabelsUp);
982  }
983  if (opt.Contains("d")) {
984  axis->SetBit(TAxis::kLabelsDown);
987  axis->ResetBit(TAxis::kLabelsUp);
988  }
989  Int_t sort = -1;
990  if (opt.Contains("a")) sort = 0;
991  if (opt.Contains(">")) sort = 1;
992  if (opt.Contains("<")) sort = 2;
993  if (sort < 0) return;
994 
995  Int_t nx = fXaxis.GetNbins()+2;
996  Int_t ny = fYaxis.GetNbins()+2;
997  Int_t n = TMath::Min(axis->GetNbins(), labels->GetSize());
998  Int_t *a = new Int_t[n+2];
999  Int_t i,j,k,bin;
1000  Double_t *sumw = new Double_t[nx*ny];
1001  Double_t *errors = new Double_t[nx*ny];
1002  Double_t *ent = new Double_t[nx*ny];
1003  THashList *labold = new THashList(labels->GetSize(),1);
1004  TIter nextold(labels);
1005  TObject *obj;
1006  while ((obj=nextold())) {
1007  labold->Add(obj);
1008  }
1009  labels->Clear();
1010  if (sort > 0) {
1011  //---sort by values of bins
1012  Double_t *pcont = new Double_t[n+2];
1013  for (i=0;i<=n;i++) pcont[i] = 0;
1014  for (i=1;i<nx;i++) {
1015  for (j=1;j<ny;j++) {
1016  bin = i+nx*j;
1017  sumw[bin] = fArray[bin];
1018  errors[bin] = fSumw2.fArray[bin];
1019  ent[bin] = fBinEntries.fArray[bin];
1020  if (axis == GetXaxis()) k = i;
1021  else k = j;
1022  if (fBinEntries.fArray[bin] != 0) pcont[k-1] += fArray[bin]/fBinEntries.fArray[bin];
1023  }
1024  }
1025  if (sort ==1) TMath::Sort(n,pcont,a,kTRUE); //sort by decreasing values
1026  else TMath::Sort(n,pcont,a,kFALSE); //sort by increasing values
1027  delete [] pcont;
1028  for (i=0;i<n;i++) {
1029  obj = labold->At(a[i]);
1030  labels->Add(obj);
1031  obj->SetUniqueID(i+1);
1032  }
1033  for (i=1;i<nx;i++) {
1034  for (j=1;j<ny;j++) {
1035  bin = i+nx*j;
1036  if (axis == GetXaxis()) {
1037  fArray[bin] = sumw[a[i-1]+1+nx*j];
1038  fSumw2.fArray[bin] = errors[a[i-1]+1+nx*j];
1039  fBinEntries.fArray[bin] = ent[a[i-1]+1+nx*j];
1040  } else {
1041  fArray[bin] = sumw[i+nx*(a[j-1]+1)];
1042  fSumw2.fArray[bin] = errors[i+nx*(a[j-1]+1)];
1043  fBinEntries.fArray[bin] = ent[i+nx*(a[j-1]+1)];
1044  }
1045  }
1046  }
1047  } else {
1048  //---alphabetic sort
1049  const UInt_t kUsed = 1<<18;
1050  TObject *objk=0;
1051  a[0] = 0;
1052  a[n+1] = n+1;
1053  for (i=1;i<=n;i++) {
1054  const char *label = "zzzzzzzzzzzz";
1055  for (j=1;j<=n;j++) {
1056  obj = labold->At(j-1);
1057  if (!obj) continue;
1058  if (obj->TestBit(kUsed)) continue;
1059  //use strcasecmp for case non-sensitive sort (may be an option)
1060  if (strcmp(label,obj->GetName()) < 0) continue;
1061  objk = obj;
1062  a[i] = j;
1063  label = obj->GetName();
1064  }
1065  if (objk) {
1066  objk->SetUniqueID(i);
1067  labels->Add(objk);
1068  objk->SetBit(kUsed);
1069  }
1070  }
1071  for (i=1;i<=n;i++) {
1072  obj = labels->At(i-1);
1073  if (!obj) continue;
1074  obj->ResetBit(kUsed);
1075  }
1076  for (i=0;i<nx;i++) {
1077  for (j=0;j<ny;j++) {
1078  bin = i+nx*j;
1079  sumw[bin] = fArray[bin];
1080  errors[bin] = fSumw2.fArray[bin];
1081  ent[bin] = fBinEntries.fArray[bin];
1082  }
1083  }
1084  for (i=0;i<nx;i++) {
1085  for (j=0;j<ny;j++) {
1086  bin = i+nx*j;
1087  if (axis == GetXaxis()) {
1088  fArray[bin] = sumw[a[i]+nx*j];
1089  fSumw2.fArray[bin] = errors[a[i]+nx*j];
1090  fBinEntries.fArray[bin] = ent[a[i]+nx*j];
1091  } else {
1092  fArray[bin] = sumw[i+nx*a[j]];
1093  fSumw2.fArray[bin] = errors[i+nx*a[j]];
1094  fBinEntries.fArray[bin] = ent[i+nx*a[j]];
1095  }
1096  }
1097  }
1098  }
1099  delete labold;
1100  if (a) delete [] a;
1101  if (sumw) delete [] sumw;
1102  if (errors) delete [] errors;
1103  if (ent) delete [] ent;
1104 }
1105 
1106 ////////////////////////////////////////////////////////////////////////////////
1107 /// Merge all histograms in the collection in this histogram.
1108 /// This function computes the min/max for the axes,
1109 /// compute a new number of bins, if necessary,
1110 /// add bin contents, errors and statistics.
1111 /// If overflows are present and limits are different the function will fail.
1112 /// The function returns the total number of entries in the result histogram
1113 /// if the merge is successful, -1 otherwise.
1114 ///
1115 /// IMPORTANT remark. The 2 axis x and y may have different number
1116 /// of bins and different limits, BUT the largest bin width must be
1117 /// a multiple of the smallest bin width and the upper limit must also
1118 /// be a multiple of the bin width.
1119 
1121 {
1122  return TProfileHelper::Merge(this, li);
1123 }
1124 
1125 ////////////////////////////////////////////////////////////////////////////////
1126 /// Performs the operation: this = this*c1*f1
1127 
1129 {
1130  Error("Multiply","Function not implemented for TProfile2D");
1131  return kFALSE;
1132 }
1133 
1134 ////////////////////////////////////////////////////////////////////////////////
1135 /// Multiply this profile2D by h1.
1136 ///
1137 /// `this = this*h1`
1138 
1140 {
1141  Error("Multiply","Multiplication of profile2D histograms not implemented");
1142  return kFALSE;
1143 }
1144 
1145 ////////////////////////////////////////////////////////////////////////////////
1146 /// Replace contents of this profile2D by multiplication of h1 by h2.
1147 ///
1148 /// `this = (c1*h1)*(c2*h2)`
1149 
1151 {
1152  Error("Multiply","Multiplication of profile2D histograms not implemented");
1153  return kFALSE;
1154 }
1155 
1156 ////////////////////////////////////////////////////////////////////////////////
1157 /// Project this profile2D into a 2-D histogram along X,Y.
1158 ///
1159 /// The projection is always of the type TH2D.
1160 ///
1161 /// - if option "E" is specified the errors of the projected histogram are computed and set
1162 /// to be equal to the errors of the profile.
1163 /// Option "E" is defined as the default one in the header file.
1164 /// - if option "" is specified the histogram errors are simply the sqrt of its content
1165 /// - if option "B" is specified, the content of bin of the returned histogram
1166 /// will be equal to the GetBinEntries(bin) of the profile,
1167 /// - if option "C=E" the bin contents of the projection are set to the
1168 /// bin errors of the profile
1169 /// - if option "W" is specified the bin content of the projected histogram is set to the
1170 /// product of the bin content of the profile and the entries.
1171 /// With this option the returned histogram will be equivalent to the one obtained by
1172 /// filling directly a TH2D using the 3-rd value as a weight.
1173 /// This option makes sense only for profile filled with all weights =1.
1174 /// When the profile is weighted (filled with weights different than 1) the
1175 /// bin error of the projected histogram (obtained using this option "W") cannot be
1176 /// correctly computed from the information stored in the profile. In that case the
1177 /// obtained histogram contains as bin error square the weighted sum of the square of the
1178 /// profiled observable (TProfile2D::fSumw2[bin] )
1179 
1180 TH2D *TProfile2D::ProjectionXY(const char *name, Option_t *option) const
1181 {
1182 
1183  TString opt = option;
1184  opt.ToLower();
1185 
1186  // Create the projection histogram
1187  // name of projected histogram is by default name of original histogram + _pxy
1188  TString pname(name);
1189  if (pname.IsNull() || pname == "_pxy")
1190  pname = TString(GetName() ) + TString("_pxy");
1191 
1192 
1193  Int_t nx = fXaxis.GetNbins();
1194  Int_t ny = fYaxis.GetNbins();
1195  const TArrayD *xbins = fXaxis.GetXbins();
1196  const TArrayD *ybins = fYaxis.GetXbins();
1197  TH2D * h1 = 0;
1198  if (xbins->fN == 0 && ybins->fN == 0) {
1199  h1 = new TH2D(pname,GetTitle(),nx,fXaxis.GetXmin(),fXaxis.GetXmax(),ny,fYaxis.GetXmin(),fYaxis.GetXmax());
1200  } else if (xbins->fN == 0) {
1201  h1 = new TH2D(pname,GetTitle(),nx,fXaxis.GetXmin(),fXaxis.GetXmax(),ny, ybins->GetArray() );
1202  } else if (ybins->fN == 0) {
1203  h1 = new TH2D(pname,GetTitle(),nx,xbins->GetArray(),ny,fYaxis.GetXmin(),fYaxis.GetXmax());
1204  } else {
1205  h1 = new TH2D(pname,GetTitle(),nx,xbins->GetArray(),ny,ybins->GetArray() );
1206  }
1207  Bool_t computeErrors = kFALSE;
1208  Bool_t cequalErrors = kFALSE;
1209  Bool_t binEntries = kFALSE;
1210  Bool_t binWeight = kFALSE;
1211  if (opt.Contains("b")) binEntries = kTRUE;
1212  if (opt.Contains("e")) computeErrors = kTRUE;
1213  if (opt.Contains("w")) binWeight = kTRUE;
1214  if (opt.Contains("c=e")) {cequalErrors = kTRUE; computeErrors=kFALSE;}
1215  if (computeErrors || binWeight || (binEntries && fBinSumw2.fN) ) h1->Sumw2();
1216 
1217  // Fill the projected histogram
1218  Int_t bin,binx, biny;
1219  Double_t cont;
1220  for (binx =0;binx<=nx+1;binx++) {
1221  for (biny =0;biny<=ny+1;biny++) {
1222  bin = GetBin(binx,biny);
1223 
1224  if (binEntries) cont = GetBinEntries(bin);
1225  else if (cequalErrors) cont = GetBinError(bin);
1226  else if (binWeight) cont = GetBinContent(bin) * GetBinEntries(bin);
1227  else cont = GetBinContent(bin); // default case
1228 
1229  h1->SetBinContent(bin ,cont);
1230 
1231  // if option E projected histogram errors are same as profile
1232  if (computeErrors ) h1->SetBinError(bin , GetBinError(bin) );
1233  // in case of option W bin error is deduced from bin sum of z**2 values of profile
1234  // this is correct only if the profile is unweighted
1235  if (binWeight) h1->GetSumw2()->fArray[bin] = fSumw2.fArray[bin];
1236  // in case of bin entries and profile is weighted, we need to set also the bin error
1237  if (binEntries && fBinSumw2.fN ) {
1238  R__ASSERT( h1->GetSumw2() );
1239  h1->GetSumw2()->fArray[bin] = fBinSumw2.fArray[bin];
1240  }
1241  }
1242  }
1243  h1->SetEntries(fEntries);
1244  return h1;
1245 }
1246 
1247 ////////////////////////////////////////////////////////////////////////////////
1248 /// Project a 2-D histogram into a profile histogram along X.
1249 ///
1250 /// The projection is made from the channels along the Y axis
1251 /// ranging from firstybin to lastybin included.
1252 /// The result is a 1D profile which contains the combination of all the considered bins along Y
1253 /// By default, bins 1 to ny are included
1254 /// When all bins are included, the number of entries in the projection
1255 /// is set to the number of entries of the 2-D histogram, otherwise
1256 /// the number of entries is incremented by 1 for all non empty cells.
1257 ///
1258 /// The option can also be used to specify the projected profile error type.
1259 /// Values which can be used are 's', 'i', or 'g'. See TProfile::BuildOptions for details
1260 
1261 TProfile *TProfile2D::ProfileX(const char *name, Int_t firstybin, Int_t lastybin, Option_t *option) const
1262 {
1263  return DoProfile(true, name, firstybin, lastybin, option);
1264 }
1265 
1266 ////////////////////////////////////////////////////////////////////////////////
1267 /// Project a 2-D histogram into a profile histogram along X
1268 ///
1269 /// The projection is made from the channels along the X axis
1270 /// ranging from firstybin to lastybin included.
1271 /// The result is a 1D profile which contains the combination of all the considered bins along X
1272 /// By default, bins 1 to ny are included
1273 /// When all bins are included, the number of entries in the projection
1274 /// is set to the number of entries of the 2-D histogram, otherwise
1275 /// the number of entries is incremented by 1 for all non empty cells.
1276 ///
1277 /// The option can also be used to specify the projected profile error type.
1278 /// Values which can be used are 's', 'i', or 'g'. See TProfile::BuildOptions for details
1279 
1280 TProfile *TProfile2D::ProfileY(const char *name, Int_t firstxbin, Int_t lastxbin, Option_t *option) const
1281 {
1282  return DoProfile(false, name, firstxbin, lastxbin, option);
1283 }
1284 
1285 ////////////////////////////////////////////////////////////////////////////////
1286 /// Implementation of ProfileX or ProfileY for a TProfile2D.
1287 ///
1288 /// Do correctly the combination of the bin averages when doing the projection
1289 
1290 TProfile * TProfile2D::DoProfile(bool onX, const char *name, Int_t firstbin, Int_t lastbin, Option_t *option) const {
1291  TString opt = option;
1292  opt.ToLower();
1293  bool originalRange = opt.Contains("o");
1294 
1295  TString expectedName = ( onX ? "_pfx" : "_pfy" );
1296 
1297  TString pname(name);
1298  if (pname.IsNull() || name == expectedName)
1299  pname = TString(GetName() ) + expectedName;
1300 
1301  const TAxis& outAxis = ( onX ? fXaxis : fYaxis );
1302  const TArrayD *bins = outAxis.GetXbins();
1303  Int_t firstOutBin = outAxis.GetFirst();
1304  Int_t lastOutBin = outAxis.GetLast();
1305 
1306  TProfile * p1 = 0;
1307  // case of fixed bins
1308  if (bins->fN == 0) {
1309  if (originalRange)
1310  p1 = new TProfile(pname,GetTitle(), outAxis.GetNbins(), outAxis.GetXmin(), outAxis.GetXmax(), opt );
1311  else
1312  p1 = new TProfile(pname,GetTitle(), lastOutBin-firstOutBin+1,
1313  outAxis.GetBinLowEdge(firstOutBin),outAxis.GetBinUpEdge(lastOutBin), opt);
1314  } else {
1315  // case of variable bins
1316  if (originalRange )
1317  p1 = new TProfile(pname,GetTitle(),outAxis.GetNbins(),bins->fArray,opt);
1318  else
1319  p1 = new TProfile(pname,GetTitle(),lastOutBin-firstOutBin+1,&bins->fArray[firstOutBin-1],opt);
1320 
1321  }
1322 
1323  if (fBinSumw2.fN) p1->Sumw2();
1324 
1325  // make projection in a 2D first
1326  TH2D * h2dW = ProjectionXY("h2temp-W","W");
1327  TH2D * h2dN = ProjectionXY("h2temp-N","B");
1328 
1329  h2dW->SetDirectory(0); h2dN->SetDirectory(0);
1330 
1331 
1332  TString opt1 = (originalRange) ? "o" : "";
1333  TH1D * h1W = (onX) ? h2dW->ProjectionX("h1temp-W",firstbin,lastbin,opt1) : h2dW->ProjectionY("h1temp-W",firstbin,lastbin,opt1);
1334  TH1D * h1N = (onX) ? h2dN->ProjectionX("h1temp-N",firstbin,lastbin,opt1) : h2dN->ProjectionY("h1temp-N",firstbin,lastbin,opt1);
1335  h1W->SetDirectory(0); h1N->SetDirectory(0);
1336 
1337 
1338  // fill the bin content
1339  R__ASSERT( h1W->fN == p1->fN );
1340  R__ASSERT( h1N->fN == p1->fN );
1341  R__ASSERT( h1W->GetSumw2()->fN != 0); // h1W should always be a weighted histogram since h2dW is
1342  for (int i = 0; i < p1->fN ; ++i) {
1343  p1->fArray[i] = h1W->GetBinContent(i); // array of profile is sum of all values
1344  p1->GetSumw2()->fArray[i] = h1W->GetSumw2()->fArray[i]; // array of content square of profile is weight square of the W projected histogram
1345  p1->SetBinEntries(i, h1N->GetBinContent(i) );
1346  if (fBinSumw2.fN) p1->GetBinSumw2()->fArray[i] = h1N->GetSumw2()->fArray[i]; // sum of weight squares are stored to compute errors in h1N histogram
1347  }
1348  // delete the created histograms
1349  delete h2dW;
1350  delete h2dN;
1351  delete h1W;
1352  delete h1N;
1353 
1354  // Also we need to set the entries since they have not been correctly calculated during the projection
1355  // we can only set them to the effective entries
1356  p1->SetEntries( p1->GetEffectiveEntries() );
1357 
1358  return p1;
1359 }
1360 
1361 
1362 ////////////////////////////////////////////////////////////////////////////////
1363 /// Replace current statistics with the values in array stats
1364 
1366 {
1367  fTsumw = stats[0];
1368  fTsumw2 = stats[1];
1369  fTsumwx = stats[2];
1370  fTsumwx2 = stats[3];
1371  fTsumwy = stats[4];
1372  fTsumwy2 = stats[5];
1373  fTsumwxy = stats[6];
1374  fTsumwz = stats[7];
1375  fTsumwz2 = stats[8];
1376 }
1377 
1378 ////////////////////////////////////////////////////////////////////////////////
1379 /// Reset contents of a Profile2D histogram.
1380 
1382 {
1383  TH2D::Reset(option);
1384  fBinEntries.Reset();
1385  fBinSumw2.Reset();
1386  TString opt = option;
1387  opt.ToUpper();
1388  if (opt.Contains("ICE") && !opt.Contains("S")) return;
1389  fTsumwz = fTsumwz2 = 0;
1390 }
1391 
1392 
1393 ////////////////////////////////////////////////////////////////////////////////
1394 /// Profile histogram is resized along axis such that x is in the axis range.
1395 ///
1396 /// The new axis limits are recomputed by doubling iteratively
1397 /// the current axis range until the specified value x is within the limits.
1398 /// The algorithm makes a copy of the histogram, then loops on all bins
1399 /// of the old histogram to fill the extended histogram.
1400 /// Takes into account errors (Sumw2) if any.
1401 /// The axis must be extendable before invoking this function.
1402 ///
1403 /// Ex: `h->GetXaxis()->SetCanExtend(kTRUE)`
1404 
1406 {
1407  TProfile2D* hold = TProfileHelper::ExtendAxis(this, x, axis);
1408  if ( hold ) {
1409  fTsumwz = hold->fTsumwz;
1410  fTsumwz2 = hold->fTsumwz2;
1411  delete hold;
1412  }
1413 }
1414 
1415 ////////////////////////////////////////////////////////////////////////////////
1416 /// Rebin this histogram grouping nxgroup/nygroup bins along the xaxis/yaxis together.
1417 ///
1418 /// if newname is not blank a new profile hnew is created.
1419 /// else the current histogram is modified (default)
1420 /// The parameter nxgroup/nygroup indicate how many bins along the xaxis/yaxis of this
1421 /// have to be merged into one bin of hnew
1422 /// If the original profile has errors stored (via Sumw2), the resulting
1423 /// profile has new errors correctly calculated.
1424 ///
1425 /// examples: if hpxpy is an existing TProfile2D profile with 40 x 40 bins
1426 /// ~~~ {.cpp}
1427 /// hpxpy->Rebin2D(); // merges two bins along the xaxis and yaxis in one
1428 /// // Carefull: previous contents of hpxpy are lost
1429 /// hpxpy->Rebin2D(3,5); // merges 3 bins along the xaxis and 5 bins along the yaxis in one
1430 /// // Carefull: previous contents of hpxpy are lost
1431 /// hpxpy->RebinX(5); //merges five bins along the xaxis in one in hpxpy
1432 /// TProfile2D *hnew = hpxpy->RebinY(5,"hnew"); // creates a new profile hnew
1433 /// // merging 5 bins of hpxpy along the yaxis in one bin
1434 /// ~~~
1435 ///
1436 /// NOTE : If nxgroup/nygroup is not an exact divider of the number of bins,
1437 /// along the xaxis/yaxis the top limit(s) of the rebinned profile
1438 /// is changed to the upper edge of the xbin=newxbins*nxgroup resp.
1439 /// ybin=newybins*nygroup and the remaining bins are added to
1440 /// the overflow bin.
1441 /// Statistics will be recomputed from the new bin contents.
1442 
1443 TProfile2D * TProfile2D::Rebin2D(Int_t nxgroup ,Int_t nygroup,const char * newname ) {
1444  //something to do?
1445  if((nxgroup != 1) || (nygroup != 1)){
1446  Int_t nxbins = fXaxis.GetNbins();
1447  Int_t nybins = fYaxis.GetNbins();
1452  if ((nxgroup <= 0) || (nxgroup > nxbins)) {
1453  Error("Rebin", "Illegal value of nxgroup=%d",nxgroup);
1454  return 0;
1455  }
1456  if ((nygroup <= 0) || (nygroup > nybins)) {
1457  Error("Rebin", "Illegal value of nygroup=%d",nygroup);
1458  return 0;
1459  }
1460 
1461  Int_t newxbins = nxbins/nxgroup;
1462  Int_t newybins = nybins/nygroup;
1463 
1464  //warning if bins are added to the overflow bin
1465  if(newxbins*nxgroup != nxbins) {
1466  Warning("Rebin", "nxgroup=%d should be an exact divider of nxbins=%d",nxgroup,nxbins);
1467  }
1468  if(newybins*nygroup != nybins) {
1469  Warning("Rebin", "nygroup=%d should be an exact divider of nybins=%d",nygroup,nybins);
1470  }
1471 
1472  //save old bin contents in new arrays
1473  Double_t *oldBins = new Double_t[(nxbins+2)*(nybins+2)];
1474  Double_t *oldCount = new Double_t[(nxbins+2)*(nybins+2)];
1475  Double_t *oldErrors = new Double_t[(nxbins+2)*(nybins+2)];
1476  Double_t *oldBinw2 = (fBinSumw2.fN ? new Double_t[(nxbins+2)*(nybins+2)] : 0 );
1477  Double_t *cu1 = GetW();
1478  Double_t *er1 = GetW2();
1479  Double_t *en1 = GetB();
1480  Double_t *ew1 = GetB2();
1481  for(Int_t ibin=0; ibin < (nxbins+2)*(nybins+2); ibin++){
1482  oldBins[ibin] = cu1[ibin];
1483  oldCount[ibin] = en1[ibin];
1484  oldErrors[ibin] = er1[ibin];
1485  if (ew1 && fBinSumw2.fN) oldBinw2[ibin] = ew1[ibin];
1486  }
1487 
1488  // create a clone of the old profile if newname is specified
1489  TProfile2D *hnew = this;
1490  if(newname && strlen(newname) > 0) {
1491  hnew = (TProfile2D*)Clone(newname);
1492  }
1493 
1494  // in case of nxgroup/nygroup not an exact divider of nxbins/nybins,
1495  // top limit is changed (see NOTE in method comment)
1496  if(newxbins*nxgroup != nxbins) {
1497  xmax = fXaxis.GetBinUpEdge(newxbins*nxgroup);
1498  hnew->fTsumw = 0; //stats must be reset because top bins will be moved to overflow bin
1499  }
1500  if(newybins*nygroup != nybins) {
1501  ymax = fYaxis.GetBinUpEdge(newybins*nygroup);
1502  hnew->fTsumw = 0; //stats must be reset because top bins will be moved to overflow bin
1503  }
1504 
1505  //rebin the axis
1506  if((fXaxis.GetXbins()->GetSize() > 0) || (fYaxis.GetXbins()->GetSize() > 0)){
1507  Double_t* xbins = new Double_t[newxbins+1];
1508  Double_t* ybins = new Double_t[newybins+1];
1509  for(Int_t i=0; i < newxbins+1; i++)
1510  xbins[i] = fXaxis.GetBinLowEdge(1+i*nxgroup);
1511  for(Int_t j=0; j < newybins+1; j++)
1512  ybins[j] = fYaxis.GetBinLowEdge(1+j*nygroup);
1513  hnew->SetBins(newxbins,xbins,newybins,ybins);
1514  delete [] xbins;
1515  delete [] ybins;
1516  }
1517  //fixed bin size
1518  else{
1519  hnew->SetBins(newxbins,xmin,xmax,newybins,ymin,ymax);
1520  }
1521 
1522  //merge bins
1523  Double_t *cu2 = hnew->GetW();
1524  Double_t *er2 = hnew->GetW2();
1525  Double_t *en2 = hnew->GetB();
1526  Double_t *ew2 = hnew->GetB2();
1527  Double_t binContent, binCount, binError, binSumw2;
1528  //connection between x and y bin number and linear global bin number:
1529  //global bin = xbin + (nxbins+2) * ybin
1530  Int_t oldxbin = 1;
1531  Int_t oldybin = 1;
1532  //global bin number
1533  Int_t bin;
1534  for(Int_t xbin = 1; xbin <= newxbins; xbin++){
1535  oldybin = 1;
1536  for(Int_t ybin = 1; ybin <= newybins; ybin++){
1537  binContent = 0;
1538  binCount = 0;
1539  binError = 0;
1540  binSumw2 = 0;
1541  for(Int_t i=0; i < nxgroup; i++){
1542  if(oldxbin + i > nxbins) break;
1543  for(Int_t j=0; j < nygroup; j++){
1544  if(oldybin + j > nybins) break;
1545  bin = oldxbin + i + (nxbins+2)*(oldybin+j);
1546  binContent += oldBins[bin];
1547  binCount += oldCount[bin];
1548  binError += oldErrors[bin];
1549  if(fBinSumw2.fN) binSumw2 += oldBinw2[bin];
1550  }
1551  }
1552  bin = xbin + (newxbins + 2)*ybin;
1553  cu2[bin] = binContent;
1554  er2[bin] = binError;
1555  en2[bin] = binCount;
1556  if(fBinSumw2.fN) ew2[bin] = binSumw2;
1557  oldybin += nygroup;
1558  }
1559  oldxbin += nxgroup;
1560  }
1561 
1562  //copy the underflow bin in x and y (0,0)
1563  cu2[0] = oldBins[0];
1564  er2[0] = oldErrors[0];
1565  en2[0] = oldCount[0];
1566  if(fBinSumw2.fN) ew2[0] = oldBinw2[0];
1567  //calculate overflow bin in x and y (newxbins+1,newybins+1)
1568  //therefore the oldxbin and oldybin from above are needed!
1569  binContent = 0;
1570  binCount = 0;
1571  binError = 0;
1572  binSumw2 = 0;
1573  for(Int_t i=oldxbin; i <= nxbins+1; i++){
1574  for(Int_t j=oldybin; j <= nybins+1; j++){
1575  //global bin number
1576  bin = i + (nxbins+2)*j;
1577  binContent += oldBins[bin];
1578  binCount += oldCount[bin];
1579  binError += oldErrors[bin];
1580  if(fBinSumw2.fN) binSumw2 += oldBinw2[bin];
1581  }
1582  }
1583  bin = (newxbins+2)*(newybins+2)-1;
1584  cu2[bin] = binContent;
1585  er2[bin] = binError;
1586  en2[bin] = binCount;
1587  if(fBinSumw2.fN) ew2[bin] = binSumw2;
1588  //calculate overflow bin in x and underflow bin in y (newxbins+1,0)
1589  binContent = 0;
1590  binCount = 0;
1591  binError = 0;
1592  binSumw2 = 0;
1593  for(Int_t i=oldxbin; i <= nxbins+1; i++){
1594  bin = i;
1595  binContent += oldBins[bin];
1596  binCount += oldCount[bin];
1597  binError += oldErrors[bin];
1598  if(fBinSumw2.fN) binSumw2 += oldBinw2[bin];
1599  }
1600  bin = newxbins + 1;
1601  cu2[bin] = binContent;
1602  er2[bin] = binError;
1603  en2[bin] = binCount;
1604  if(fBinSumw2.fN) ew2[bin] = binSumw2;
1605  //calculate underflow bin in x and overflow bin in y (0,newybins+1)
1606  binContent = 0;
1607  binCount = 0;
1608  binError = 0;
1609  binSumw2 = 0;
1610  for(Int_t i=oldybin; i <= nybins+1; i++){
1611  bin = i*(nxbins + 2);
1612  binContent += oldBins[bin];
1613  binCount += oldCount[bin];
1614  binError += oldErrors[bin];
1615  if(fBinSumw2.fN) binSumw2 += oldBinw2[bin];
1616  }
1617  bin = (newxbins + 2)*(newybins + 1);
1618  cu2[bin] = binContent;
1619  er2[bin] = binError;
1620  en2[bin] = binCount;
1621  if(fBinSumw2.fN) ew2[bin] = binSumw2;
1622  //calculate under/overflow contents in y for the new x bins
1623  Double_t binContentuf, binCountuf, binErroruf, binSumw2uf;
1624  Double_t binContentof, binCountof, binErrorof, binSumw2of;
1625  Int_t ufbin, ofbin;
1626  Int_t oldxbin2 = 1;
1627  for(Int_t xbin = 1; xbin <= newxbins; xbin++){
1628  binContentuf = 0;
1629  binCountuf = 0;
1630  binErroruf = 0;
1631  binSumw2uf = 0;
1632  binContentof = 0;
1633  binCountof = 0;
1634  binErrorof = 0;
1635  binSumw2of = 0;
1636  for(Int_t i = 0; i < nxgroup; i++){
1637  //index of under/overflow bin for y in old binning
1638  ufbin = (oldxbin2 + i);
1639  binContentuf += oldBins[ufbin];
1640  binCountuf += oldCount[ufbin];
1641  binErroruf += oldErrors[ufbin];
1642  if(fBinSumw2.fN) binSumw2uf += oldBinw2[ufbin];
1643  for(Int_t j = oldybin; j <= nybins+1; j++)
1644  {
1645  ofbin = ufbin + j*(nxbins + 2);
1646  binContentof += oldBins[ofbin];
1647  binCountof += oldCount[ofbin];
1648  binErrorof += oldErrors[ofbin];
1649  if(fBinSumw2.fN) binSumw2of += oldBinw2[ofbin];
1650  }
1651  }
1652  //index of under/overflow bin for y in new binning
1653  ufbin = xbin;
1654  ofbin = ufbin + (newybins + 1)*(newxbins + 2);
1655  cu2[ufbin] = binContentuf;
1656  er2[ufbin] = binErroruf;
1657  en2[ufbin] = binCountuf;
1658  if(fBinSumw2.fN) ew2[ufbin] = binSumw2uf;
1659  cu2[ofbin] = binContentof;
1660  er2[ofbin] = binErrorof;
1661  en2[ofbin] = binCountof;
1662  if(fBinSumw2.fN) ew2[ofbin] = binSumw2of;
1663 
1664  oldxbin2 += nxgroup;
1665  }
1666  //calculate under/overflow contents in x for the new y bins
1667  Int_t oldybin2 = 1;
1668  for(Int_t ybin = 1; ybin <= newybins; ybin++){
1669  binContentuf = 0;
1670  binCountuf = 0;
1671  binErroruf = 0;
1672  binSumw2uf = 0;
1673  binContentof = 0;
1674  binCountof = 0;
1675  binErrorof = 0;
1676  binSumw2of = 0;
1677  for(Int_t i = 0; i < nygroup; i++){
1678  //index of under/overflow bin for x in old binning
1679  ufbin = (oldybin2 + i)*(nxbins+2);
1680  binContentuf += oldBins[ufbin];
1681  binCountuf += oldCount[ufbin];
1682  binErroruf += oldErrors[ufbin];
1683  if(fBinSumw2.fN) binSumw2uf += oldBinw2[ufbin];
1684  for(Int_t j = oldxbin; j <= nxbins+1; j++)
1685  {
1686  ofbin = j + ufbin;
1687  binContentof += oldBins[ofbin];
1688  binCountof += oldCount[ofbin];
1689  binErrorof += oldErrors[ofbin];
1690  if(fBinSumw2.fN) binSumw2of += oldBinw2[ofbin];
1691  }
1692  }
1693  //index of under/overflow bin for x in new binning
1694  ufbin = ybin * (newxbins + 2);
1695  ofbin = newxbins + 1 + ufbin;
1696  cu2[ufbin] = binContentuf;
1697  er2[ufbin] = binErroruf;
1698  en2[ufbin] = binCountuf;
1699  if(fBinSumw2.fN) ew2[ufbin] = binSumw2uf;
1700  cu2[ofbin] = binContentof;
1701  er2[ofbin] = binErrorof;
1702  en2[ofbin] = binCountof;
1703  if(fBinSumw2.fN) ew2[ofbin] = binSumw2of;
1704 
1705  oldybin2 += nygroup;
1706  }
1707 
1708  delete [] oldBins;
1709  delete [] oldCount;
1710  delete [] oldErrors;
1711  if (oldBinw2) delete [] oldBinw2;
1712 
1713  return hnew;
1714  }
1715  //nxgroup == nygroup == 1
1716  else{
1717  if((newname) && (strlen(newname) > 0))
1718  return (TProfile2D*)Clone(newname);
1719  else
1720  return this;
1721  }
1722 }
1723 
1724 ////////////////////////////////////////////////////////////////////////////////
1725 /// Rebin only the X axis.
1726 /// see Rebin2D
1727 
1728 TProfile2D * TProfile2D::RebinX(Int_t ngroup,const char * newname ) {
1729  return Rebin2D(ngroup,1,newname);
1730 }
1731 
1732 ////////////////////////////////////////////////////////////////////////////////
1733 /// Rebin only the Y axis.
1734 /// see Rebin2D
1735 
1736 TProfile2D * TProfile2D::RebinY(Int_t ngroup,const char * newname ) {
1737  return Rebin2D(1,ngroup,newname);
1738 }
1739 
1740 ////////////////////////////////////////////////////////////////////////////////
1741 /// Save primitive as a C++ statement(s) on output stream out.
1742 ///
1743 /// Note the following restrictions in the code generated:
1744 /// - variable bin size not implemented
1745 /// - SetErrorOption not implemented
1746 
1747 void TProfile2D::SavePrimitive(std::ostream &out, Option_t *option /*= ""*/)
1748 {
1749  char quote = '"';
1750  out <<" "<<std::endl;
1751  out <<" "<<ClassName()<<" *";
1752 
1753  out << GetName() << " = new " << ClassName() << "(" << quote
1754  << GetName() << quote << "," << quote<< GetTitle() << quote
1755  << "," << GetXaxis()->GetNbins();
1756  out << "," << GetXaxis()->GetXmin()
1757  << "," << GetXaxis()->GetXmax();
1758  out << "," << GetYaxis()->GetNbins();
1759  out << "," << GetYaxis()->GetXmin()
1760  << "," << GetYaxis()->GetXmax();
1761  out << "," << fZmin
1762  << "," << fZmax;
1763  out << ");" << std::endl;
1764 
1765 
1766  // save bin entries
1767  Int_t bin;
1768  for (bin=0;bin<fNcells;bin++) {
1769  Double_t bi = GetBinEntries(bin);
1770  if (bi) {
1771  out<<" "<<GetName()<<"->SetBinEntries("<<bin<<","<<bi<<");"<<std::endl;
1772  }
1773  }
1774  //save bin contents
1775  for (bin=0;bin<fNcells;bin++) {
1776  Double_t bc = fArray[bin];
1777  if (bc) {
1778  out<<" "<<GetName()<<"->SetBinContent("<<bin<<","<<bc<<");"<<std::endl;
1779  }
1780  }
1781  // save bin errors
1782  if (fSumw2.fN) {
1783  for (bin=0;bin<fNcells;bin++) {
1784  Double_t be = TMath::Sqrt(fSumw2.fArray[bin]);
1785  if (be) {
1786  out<<" "<<GetName()<<"->SetBinError("<<bin<<","<<be<<");"<<std::endl;
1787  }
1788  }
1789  }
1790 
1791  TH1::SavePrimitiveHelp(out, GetName(), option);
1792 }
1793 
1794 ////////////////////////////////////////////////////////////////////////////////
1795 /// Multiply this profile2D by a constant c1.
1796 ///
1797 /// `this = c1*this
1798 ///
1799 /// This function uses the services of TProfile2D::Add
1800 
1802 {
1803  TProfileHelper::Scale(this, c1, option);
1804 }
1805 
1806 ////////////////////////////////////////////////////////////////////////////////
1807 /// Set the number of entries in bin.
1808 
1810 {
1811  TProfileHelper::SetBinEntries(this, bin, w);
1812 }
1813 
1814 ////////////////////////////////////////////////////////////////////////////////
1815 /// Redefine x and y axis parameters.
1816 
1818 {
1819  TH1::SetBins(nx,xmin, xmax,ny, ymin,ymax);
1822 }
1823 
1824 ////////////////////////////////////////////////////////////////////////////////
1825 /// Redefine x and y axis parameters for variable bin sizes.
1826 
1827 void TProfile2D::SetBins(Int_t nx, const Double_t *xbins, Int_t ny, const Double_t *ybins)
1828 {
1829  TH1::SetBins(nx,xbins,ny,ybins);
1832 }
1833 
1834 ////////////////////////////////////////////////////////////////////////////////
1835 /// Set total number of bins including under/overflow.
1836 /// Reallocate bin contents array
1837 
1839 {
1842 }
1843 
1844 ////////////////////////////////////////////////////////////////////////////////
1845 /// Set the buffer size in units of 8 bytes (double).
1846 
1848 {
1849  if (fBuffer) {
1850  BufferEmpty();
1851  delete [] fBuffer;
1852  fBuffer = 0;
1853  }
1854  if (buffersize <= 0) {
1855  fBufferSize = 0;
1856  return;
1857  }
1858  if (buffersize < 100) buffersize = 100;
1859  fBufferSize = 1 + 4*buffersize;
1860  fBuffer = new Double_t[fBufferSize];
1861  memset(fBuffer,0,sizeof(Double_t)*fBufferSize);
1862 }
1863 
1864 ////////////////////////////////////////////////////////////////////////////////
1865 /// Set option to compute profile2D errors.
1866 ///
1867 /// The computation of the bin errors is based on the parameter option:
1868 /// - ' ' (Default) The bin errors are the standard error on the mean of the bin profiled values (Z),
1869 /// i.e. the standard error of the bin contents.
1870 /// Note that if TProfile::Approximate() is called, an approximation is used when
1871 /// the spread in Z is 0 and the number of bin entries is > 0
1872 /// - 's' The bin errors are the standard deviations of the Z bin values
1873 /// Note that if TProfile::Approximate() is called, an approximation is used when
1874 /// the spread in Z is 0 and the number of bin entries is > 0
1875 /// - 'i' Errors are as in default case (standard errors of the bin contents)
1876 /// The only difference is for the case when the spread in Z is zero.
1877 /// In this case for N > 0 the error is 1./SQRT(12.*N)
1878 /// - 'g' Errors are 1./SQRT(W) for W not equal to 0 and 0 for W = 0.
1879 /// W is the sum in the bin of the weights of the profile.
1880 /// This option is for combining measurements z +/- dz,
1881 /// and the profile is filled with values y and weights z = 1/dz**2
1882 ///
1883 /// See TProfile::BuildOptions for a detailed explanation of all options
1884 
1886 {
1887  TProfileHelper::SetErrorOption(this, option);
1888 }
1889 
1890 ////////////////////////////////////////////////////////////////////////////////
1891 /// Stream an object of class TProfile2D.
1892 
1893 void TProfile2D::Streamer(TBuffer &R__b)
1894 {
1895  if (R__b.IsReading()) {
1896  UInt_t R__s, R__c;
1897  Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
1898  if (R__v > 2) {
1899  R__b.ReadClassBuffer(TProfile2D::Class(), this, R__v, R__s, R__c);
1900  return;
1901  }
1902  //====process old versions before automatic schema evolution
1903  TH2D::Streamer(R__b);
1904  fBinEntries.Streamer(R__b);
1905  Int_t errorMode;
1906  R__b >> errorMode;
1907  fErrorMode = (EErrorType)errorMode;
1908  if (R__v < 2) {
1909  Float_t zmin,zmax;
1910  R__b >> zmin; fZmin = zmin;
1911  R__b >> zmax; fZmax = zmax;
1912  } else {
1913  R__b >> fZmin;
1914  R__b >> fZmax;
1915  }
1916  R__b.CheckByteCount(R__s, R__c, TProfile2D::IsA());
1917  //====end of old versions
1918 
1919  } else {
1920  R__b.WriteClassBuffer(TProfile2D::Class(),this);
1921  }
1922 }
1923 
1924 ////////////////////////////////////////////////////////////////////////////////
1925 /// Create/Delete structure to store sum of squares of weights per bin.
1926 ///
1927 /// This is needed to compute the correct statistical quantities
1928 /// of a profile filled with weights
1929 ///
1930 /// This function is automatically called when the histogram is created
1931 /// if the static function TH1::SetDefaultSumw2 has been called before.
1932 /// If flag is false the structure is deleted
1933 
1935 {
1936  TProfileHelper::Sumw2(this, flag);
1937 }
const int nx
Definition: kalman.C:16
TArrayD()
Default TArrayD ctor.
Definition: TArrayD.cxx:26
virtual void SetBinsLength(Int_t n=-1)
Set total number of bins including under/overflow Reallocate bin contents array.
Definition: TH2.cxx:3818
virtual Long64_t Merge(TCollection *list)
Merge all histograms in the collection in this histogram.
virtual Double_t GetBinEntries(Int_t bin) const
Return bin entries of a Profile2D histogram.
Definition: TProfile2D.cxx:791
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
static void SetBinEntries(T *p, Int_t bin, Double_t w)
Bool_t IsReading() const
Definition: TBuffer.h:81
static void Scale(T *p, Double_t c1, Option_t *option)
virtual Double_t GetEffectiveEntries() const
Number of effective entries of the histogram.
Definition: TH1.cxx:4079
float xmin
Definition: THbookFile.cxx:93
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
Int_t GetFirst() const
Return first bin on the axis i.e.
Definition: TAxis.cxx:444
long long Long64_t
Definition: RtypesCore.h:69
virtual void Sumw2(Bool_t flag=kTRUE)
Create/Delete structure to store sum of squares of weights per bin.
short Version_t
Definition: RtypesCore.h:61
static Double_t GetBinError(T *p, Int_t bin)
virtual Double_t GetBinError(Int_t bin) const
Return bin error of a Profile2D histogram.
Definition: TProfile2D.cxx:830
float Float_t
Definition: RtypesCore.h:53
virtual void SetDirectory(TDirectory *dir)
By default when an histogram is created, it is added to the list of histogram objects in the current ...
Definition: TH1.cxx:8053
const char Option_t
Definition: RtypesCore.h:62
virtual Double_t GetBinLowEdge(Int_t bin) const
Return low edge of bin.
Definition: TAxis.cxx:504
return c1
Definition: legend1.C:41
void Reset()
Definition: TArrayD.h:47
virtual void ExtendAxis(Double_t x, TAxis *axis)
Profile histogram is resized along axis such that x is in the axis range.
float ymin
Definition: THbookFile.cxx:93
TProfile * ProfileY(const char *name="_pfy", Int_t firstxbin=0, Int_t lastxbin=-1, Option_t *option="") const
Project a 2-D histogram into a profile histogram along X.
Option_t * GetErrorOption() const
Return option to compute profile2D errors.
Definition: TProfile2D.cxx:838
TAxis fYaxis
Y axis descriptor.
Definition: TH1.h:81
static void SetErrorOption(T *p, Option_t *opt)
TH1D * ProjectionY(const char *name="_py", Int_t firstxbin=0, Int_t lastxbin=-1, Option_t *option="") const
Project a 2-D histogram into a 1-D histogram along Y.
Definition: TH2.cxx:2328
const Double_t * GetArray() const
Definition: TArrayD.h:43
Bool_t TestBit(UInt_t f) const
Definition: TObject.h:159
virtual void SetBins(Int_t nx, Double_t xmin, Double_t xmax)
Redefine x axis parameters.
Definition: TH1.cxx:7883
static void LabelsInflate(T *p, Option_t *)
static Bool_t fgStatOverflows
!flag to use under/overflows in statistics
Definition: TH1.h:106
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH1.cxx:4639
void ToUpper()
Change string to upper case.
Definition: TString.cxx:1112
Buffer base class used for serializing objects.
Definition: TBuffer.h:40
TH1D * ProjectionX(const char *name="_px", Int_t firstybin=0, Int_t lastybin=-1, Option_t *option="") const
Project a 2-D histogram into a 1-D histogram along X.
Definition: TH2.cxx:2288
#define R__ASSERT(e)
Definition: TError.h:96
static THLimitsFinder * GetLimitsFinder()
Return pointer to the current finder.
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)=0
virtual TProfile2D * RebinY(Int_t ngroup=2, const char *newname="")
Rebin only the Y axis.
Basic string class.
Definition: TString.h:129
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:168
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1099
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
Histogram is forced to be not weighted even when the histogram is filled with weighted different than...
Definition: TH1.h:155
virtual TProfile2D * Rebin2D(Int_t nxgroup=2, Int_t nygroup=2, const char *newname="")
Rebin this histogram grouping nxgroup/nygroup bins along the xaxis/yaxis together.
TArrayD fSumw2
Array of sum of squares of weights.
Definition: TH1.h:94
virtual ~TProfile2D()
Default destructor for Profile2D histograms.
Definition: TProfile2D.cxx:83
Profile Histogram.
Definition: TProfile.h:32
virtual Int_t FindGoodLimits(TH1 *h, Double_t xmin, Double_t xmax)
compute the best axis limits for the X axis.
Double_t * GetB()
Definition: TProfile2D.h:62
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save primitive as a C++ statement(s) on output stream out.
Short_t Abs(Short_t d)
Definition: TMathBase.h:108
virtual Bool_t Multiply(TF1 *h1, Double_t c1=1)
Performs the operation: this = this*c1*f1.
EErrorType
Definition: TProfile.h:28
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:687
Double_t fTsumwx2
Total Sum of weight*X*X.
Definition: TH1.h:89
virtual Bool_t CanExtendAllAxes() const
Returns true if all axes are extendable.
Definition: TH1.cxx:5961
virtual Double_t GetBinUpEdge(Int_t bin) const
Return up edge of bin.
Definition: TAxis.cxx:514
static void BuildArray(T *p)
static Bool_t fgApproximate
Definition: TProfile2D.h:41
Int_t Fill(const Double_t *v)
Definition: TProfile2D.h:50
Double_t GetXmin() const
Definition: TAxis.h:133
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition: TObject.cxx:135
Double_t x[n]
Definition: legend1.C:17
virtual Bool_t Divide(TF1 *h1, Double_t c1=1)
Performs the operation: this = this/(c1*f1) .
Definition: TProfile2D.cxx:393
void Class()
Definition: Class.C:29
THashList implements a hybrid collection class consisting of a hash table and a list to store TObject...
Definition: THashList.h:34
void BuildOptions(Double_t zmin, Double_t zmax, Option_t *option)
Set Profile2D histogram structure and options.
Definition: TProfile2D.cxx:174
const int ny
Definition: kalman.C:17
Double_t * GetB2()
Definition: TProfile2D.h:63
TArrayD fBinEntries
Definition: TProfile2D.h:33
virtual Double_t GetBinContent(Int_t bin) const
Return bin content of a Profile2D histogram.
Definition: TProfile2D.cxx:778
THashList * GetLabels() const
Definition: TAxis.h:117
virtual void SetErrorOption(Option_t *option="")
Set option to compute profile2D errors.
virtual TArrayD * GetBinSumw2()
Definition: TProfile.h:109
static double p2(double t, double a, double b, double c)
static void LabelsDeflate(T *p, Option_t *)
static void Approximate(Bool_t approx=kTRUE)
Static function, set the fgApproximate flag.
Definition: TProfile2D.cxx:255
virtual Double_t GetBinCenter(Int_t bin) const
Return center of bin.
Definition: TAxis.cxx:464
void Sort(Index n, const Element *a, Index *index, Bool_t down=kTRUE)
Definition: TMath.h:1151
Double_t * fArray
Definition: TArrayD.h:30
static T * ExtendAxis(T *p, Double_t x, TAxis *axis)
Double_t fTsumwy2
Definition: TH2.h:35
virtual void LabelsDeflate(Option_t *axis="X")
Reduce the number of bins for this axis to the number of bins having a label.
Definition: TProfile2D.cxx:925
virtual TArrayD * GetSumw2()
Definition: TH1.h:293
TH1F * h1
Definition: legend1.C:5
virtual void SetBinError(Int_t bin, Double_t error)
See convention for numbering bins in TH1::GetBin.
Definition: TH1.cxx:8311
virtual void SetBinsLength(Int_t n=-1)
Set total number of bins including under/overflow.
virtual void Scale(Double_t c1=1, Option_t *option="")
Multiply this profile2D by a constant c1.
Double_t fTsumwx
Total Sum of weight*X.
Definition: TH1.h:88
Double_t fZmin
Definition: TProfile2D.h:35
virtual void SetUniqueID(UInt_t uid)
Set the unique object id.
Definition: TObject.cxx:698
virtual void AddBinContent(Int_t bin)
Increment bin content by 1.
Definition: TH2.h:307
Int_t fN
Definition: TArray.h:38
void Clear(Option_t *option="")
Remove all objects from the list.
Definition: THashList.cxx:168
virtual Bool_t Add(TF1 *h1, Double_t c1=1, Option_t *option="")
Performs the operation: this = this + c1*f1 .
Definition: TProfile2D.cxx:199
float ymax
Definition: THbookFile.cxx:93
Int_t GetLast() const
Return last bin on the axis i.e.
Definition: TAxis.cxx:455
Class to manage histogram axis.
Definition: TAxis.h:30
TH2D * ProjectionXY(const char *name="_pxy", Option_t *option="e") const
Project this profile2D into a 2-D histogram along X,Y.
Int_t GetSize() const
Definition: TArray.h:47
Double_t fTsumwxy
Definition: TH2.h:36
Double_t fTsumwy
Definition: TH2.h:34
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:436
Collection abstract base class.
Definition: TCollection.h:42
virtual Int_t BufferEmpty(Int_t action=0)
Fill histogram with all entries in the buffer.
Definition: TProfile2D.cxx:269
static Int_t fgBufferSize
!default buffer size for automatic histograms
Definition: TH1.h:104
virtual void Copy(TObject &hnew) const
Copy a Profile2D histogram to a new profile2D histogram.
Definition: TProfile2D.cxx:364
virtual void GetStats(Double_t *stats) const
Fill the array stats from the contents of this profile.
Definition: TProfile2D.cxx:865
virtual void Copy(TObject &hnew) const
Copy.
Definition: TH2.cxx:3798
Double_t fZmax
Definition: TProfile2D.h:36
virtual void PutStats(Double_t *stats)
Replace current statistics with the values in array stats.
Double_t fEntries
Number of entries.
Definition: TH1.h:85
unsigned int UInt_t
Definition: RtypesCore.h:42
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:873
Double_t fTsumwz2
Definition: TProfile2D.h:39
virtual TProfile2D * RebinX(Int_t ngroup=2, const char *newname="")
Rebin only the X axis.
virtual void LabelsOption(Option_t *option="h", Option_t *axis="X")
Set option(s) to draw axis with labels.
Definition: TProfile2D.cxx:953
TAxis * GetYaxis()
Definition: TH1.h:301
static double p1(double t, double a, double b)
float xmax
Definition: THbookFile.cxx:93
virtual TObject * At(Int_t idx) const
Returns the object at position idx. Returns 0 if idx is out of range.
Definition: TList.cxx:315
virtual void SetBuffer(Int_t buffersize, Option_t *option="")
Set the buffer size in units of 8 bytes (double).
virtual void LabelsInflate(Option_t *axis="X")
Double the number of bins for axis.
Definition: TProfile2D.cxx:935
tomato 1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:594
TProfile2D()
Default constructor for Profile2D histograms.
Definition: TProfile2D.cxx:73
const Bool_t kFALSE
Definition: RtypesCore.h:92
Double_t * GetW()
Definition: TProfile2D.h:64
virtual Int_t FindBin(Double_t x)
Find bin number corresponding to abscissa x.
Definition: TAxis.cxx:279
static Bool_t Add(T *p, const TH1 *h1, const TH1 *h2, Double_t c1, Double_t c2=1)
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
Double_t fTsumw2
Total Sum of squares of weights.
Definition: TH1.h:87
TProfile * ProfileX(const char *name="_pfx", Int_t firstybin=0, Int_t lastybin=-1, Option_t *option="") const
Project a 2-D histogram into a profile histogram along X.
virtual Int_t BufferFill(Double_t, Double_t)
accumulate arguments in buffer.
Definition: TProfile2D.h:43
TH2D()
Constructor.
Definition: TH2.cxx:3689
return c2
Definition: legend2.C:14
#define ClassImp(name)
Definition: Rtypes.h:336
double Double_t
Definition: RtypesCore.h:55
virtual void SavePrimitiveHelp(std::ostream &out, const char *hname, Option_t *option="")
Helper function for the SavePrimitive functions from TH1 or classes derived from TH1, eg TProfile, TProfile2D.
Definition: TH1.cxx:6614
virtual Double_t GetBinEffectiveEntries(Int_t bin)
Return bin effective entries for a weighted filled Profile histogram.
Definition: TProfile2D.cxx:807
Double_t fTsumw
Total Sum of weights.
Definition: TH1.h:86
static Double_t GetBinEffectiveEntries(T *p, Int_t bin)
Double_t y[n]
Definition: legend1.C:17
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:572
The TH1 histogram class.
Definition: TH1.h:56
virtual void SetBinEntries(Int_t bin, Double_t w)
Set the number of entries in bin.
Profile2D histograms are used to display the mean value of Z and its RMS for each cell in X...
Definition: TProfile2D.h:27
Array of doubles (64 bits per element).
Definition: TArrayD.h:27
Bool_t IsNull() const
Definition: TString.h:385
Mother of all ROOT objects.
Definition: TObject.h:37
virtual void SetBinEntries(Int_t bin, Double_t w)
Set the number of entries in bin.
Definition: TProfile.cxx:1620
you should not use this method at all Int_t Int_t z
Definition: TRolke.cxx:630
virtual void Add(TObject *obj)
Definition: TList.h:77
Int_t fBufferSize
fBuffer size
Definition: TH1.h:97
1-Dim function class
Definition: TF1.h:150
Int_t IsNaN(Double_t x)
Definition: TMath.h:778
virtual void Sumw2(Bool_t flag=kTRUE)
Create structure to store sum of squares of weights.
Definition: TH1.cxx:8132
TObject * Clone(const char *newname=0) const
Make a complete copy of the underlying object.
Definition: TH1.cxx:2544
TArrayD fBinSumw2
Definition: TProfile2D.h:40
virtual void SetEntries(Double_t n)
Definition: TH1.h:363
TAxis fXaxis
X axis descriptor.
Definition: TH1.h:80
static void Sumw2(T *p, Bool_t flag)
void ResetBit(UInt_t f)
Definition: TObject.h:158
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content.
Definition: TH2.cxx:2471
virtual Int_t GetNbinsX() const
Definition: TH1.h:277
Double_t Sqrt(Double_t x)
Definition: TMath.h:591
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition: TObject.cxx:901
Int_t GetNbins() const
Definition: TAxis.h:121
virtual void Sumw2(Bool_t flag=kTRUE)
Create/delete structure to store sum of squares of weights per bin.
Definition: TProfile.cxx:1745
virtual const char * GetName() const
Returns name of object.
Definition: TObject.cxx:364
virtual Int_t GetSize() const
Definition: TCollection.h:89
virtual Int_t GetBin(Int_t binx, Int_t biny, Int_t binz=0) const
Return Global bin number corresponding to binx,y,z.
Definition: TH2.cxx:966
void SetBins(const Int_t *nbins, const Double_t *range)
Definition: TProfile2D.h:48
Double_t * fBuffer
[fBufferSize] entry buffer
Definition: TH1.h:98
EErrorType fErrorMode
Definition: TProfile2D.h:34
const Bool_t kTRUE
Definition: RtypesCore.h:91
void Set(Int_t n)
Set size of this array to n doubles.
Definition: TArrayD.cxx:105
Double_t GetXmax() const
Definition: TAxis.h:134
virtual TProfile * DoProfile(bool onX, const char *name, Int_t firstbin, Int_t lastbin, Option_t *option) const
Implementation of ProfileX or ProfileY for a TProfile2D.
Double_t fTsumwz
True when TProfile2D::Scale is called.
Definition: TProfile2D.h:38
void Copy(TArrayD &array) const
Definition: TArrayD.h:42
const Int_t n
Definition: legend1.C:16
static Long64_t Merge(T *p, TCollection *list)
Double_t * GetW2()
Definition: TProfile2D.h:65
const TArrayD * GetXbins() const
Definition: TAxis.h:130
Bool_t fScaling
Definition: TProfile2D.h:37
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:859
TAxis * GetXaxis()
Definition: TH1.h:300
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
virtual Int_t GetNbinsY() const
Definition: TH1.h:278
Int_t fNcells
number of bins(1D), cells (2D) +U/Overflows
Definition: TH1.h:79
tomato 2-D histogram with a double per channel (see TH1 documentation)}
Definition: TH2.h:290