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