Logo ROOT   6.10/09
Reference Guide
TProfile3D.cxx
Go to the documentation of this file.
1 // @(#)root/hist:$Id$
2 // Author: Rene Brun 17/05/2006
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 "TProfile3D.h"
13 #include "TProfile2D.h"
14 #include "THashList.h"
15 #include "TMath.h"
16 #include "THLimitsFinder.h"
17 #include "Riostream.h"
18 #include "TVirtualPad.h"
19 #include "TError.h"
20 #include "TClass.h"
21 
22 #include "TProfileHelper.h"
23 
25 
27 
28 /** \class TProfile3D
29  \ingroup Hist
30  Profile3D histograms are used to display the mean
31  value of T and its RMS for each cell in X,Y,Z.
32  Profile3D histograms are in many cases an
33  The inter-relation of three measured quantities X, Y, Z and T can always
34  be visualized by a four-dimensional histogram or scatter-plot;
35  its representation on the line-printer is not particularly
36  satisfactory, except for sparse data. If T is an unknown (but single-valued)
37  approximate function of X,Y,Z this function is displayed by a profile3D histogram with
38  much better precision than by a scatter-plot.
39 
40  The following formulae show the cumulated contents (capital letters) and the values
41  displayed by the printing or plotting routines (small letters) of the elements for cell I, J.
42 
43  2
44  H(I,J,K) = sum T E(I,J,K) = sum T
45  l(I,J,K) = sum l L(I,J,K) = sum l
46  h(I,J,K) = H(I,J,K)/L(I,J,K) s(I,J,K) = sqrt(E(I,J,K)/L(I,J,K)- h(I,J,K)**2)
47  e(I,J,K) = s(I,J,K)/sqrt(L(I,J,K))
48 
49  In the special case where s(I,J,K) is zero (eg, case of 1 entry only in one cell)
50  e(I,J,K) is computed from the average of the s(I,J,K) for all cells,
51  if the static function TProfile3D::Approximate has been called.
52  This simple/crude approximation was suggested in order to keep the cell
53  during a fit operation. But note that this approximation is not the default behaviour.
54 
55  Example of a profile3D histogram
56 ~~~~{.cpp}
57 {
58  auto c1 = new TCanvas("c1","Profile histogram example",200,10,700,500);
59  auto hprof3d = new TProfile3D("hprof3d","Profile of pt versus px, py and pz",40,-4,4,40,-4,4,40,0,20);
60  Double_t px, py, pz, pt;
61  TRandom3 r(0);
62  for ( Int_t i=0; i<25000; i++) {
63  r.Rannor(px,py);
64  pz = px*px + py*py;
65  pt = r.Landau(0,1);
66  hprof3d->Fill(px,py,pz,pt,1);
67  }
68  hprof3d->Draw();
69 }
70 ~~~~
71  NOTE: A TProfile3D is drawn as it was a simple TH3
72 */
73 
74 ////////////////////////////////////////////////////////////////////////////////
75 /// Default constructor for Profile3D histograms.
76 
78 {
79  fTsumwt = fTsumwt2 = 0;
80  fScaling = kFALSE;
81  BuildOptions(0,0,"");
82 }
83 
84 ////////////////////////////////////////////////////////////////////////////////
85 /// Default destructor for Profile3D histograms.
86 
88 {
89 }
90 
91 ////////////////////////////////////////////////////////////////////////////////
92 /// Normal Constructor for Profile histograms.
93 ///
94 /// The first eleven parameters are similar to TH3D::TH3D.
95 /// All values of t are accepted at filling time.
96 /// To fill a profile3D histogram, one must use TProfile3D::Fill function.
97 ///
98 /// Note that when filling the profile histogram the function Fill
99 /// checks if the variable t is between fTmin and fTmax.
100 /// If a minimum or maximum value is set for the T scale before filling,
101 /// then all values below tmin or above tmax will be discarded.
102 /// Setting the minimum or maximum value for the T scale before filling
103 /// has the same effect as calling the special TProfile3D constructor below
104 /// where tmin and tmax are specified.
105 ///
106 /// H(I,J,K) is printed as the cell contents. The errors computed are s(I,J,K) if CHOPT='S'
107 /// (spread option), or e(I,J,K) if CHOPT=' ' (error on mean).
108 ///
109 /// See TProfile3D::BuildOptions for explanation of parameters
110 ///
111 /// see other constructors below with all possible combinations of
112 /// fix and variable bin size like in TH3D.
113 
114 TProfile3D::TProfile3D(const char *name,const char *title,Int_t nx,Double_t xlow,Double_t xup,Int_t ny,Double_t ylow,Double_t yup,Int_t nz, Double_t zlow,Double_t zup,Option_t *option)
115  : TH3D(name,title,nx,xlow,xup,ny,ylow,yup,nz,zlow,zup)
116 {
117  BuildOptions(0,0,option);
118  if (xlow >= xup || ylow >= yup || zlow >= zup) SetBuffer(fgBufferSize);
119 }
120 
121 ////////////////////////////////////////////////////////////////////////////////
122 /// Create a 3-D Profile with variable bins in X , Y and Z.
123 
124 TProfile3D::TProfile3D(const char *name,const char *title,Int_t nx,const Double_t *xbins,Int_t ny,const Double_t *ybins,Int_t nz,const Double_t *zbins,Option_t *option)
125  : TH3D(name,title,nx,xbins,ny,ybins,nz,zbins)
126 {
127  BuildOptions(0,0,option);
128 }
129 
130 ////////////////////////////////////////////////////////////////////////////////
131 /// Set Profile3D histogram structure and options.
132 ///
133 /// - tmin: minimum value allowed for t
134 /// - tmax: maximum value allowed for t
135 /// if (tmin = tmax = 0) there are no limits on the allowed t values (tmin = -inf, tmax = +inf)
136 ///
137 /// - option: this is the option for the computation of the t error of the profile ( TProfile3D::GetBinError )
138 /// possible values for the options are documented in TProfile3D::SetErrorOption
139 ///
140 /// see also TProfile::BuildOptions for a detailed description
141 
143 {
144  SetErrorOption(option);
145 
146  // create extra profile data structure (bin entries/ y^2 and sum of weight square)
148 
149  fTmin = tmin;
150  fTmax = tmax;
151  fScaling = kFALSE;
152  fTsumwt = fTsumwt2 = 0;
153 }
154 
155 ////////////////////////////////////////////////////////////////////////////////
156 /// Copy constructor.
157 
159 {
160  ((TProfile3D&)profile).Copy(*this);
161 }
162 
163 ////////////////////////////////////////////////////////////////////////////////
164 /// Performs the operation: `this = this + c1*f1` .
165 
167 {
168  Error("Add","Function not implemented for TProfile3D");
169  return kFALSE;
170 }
171 
172 ////////////////////////////////////////////////////////////////////////////////
173 /// Performs the operation: `this = this + c1*h1` .
174 
176 {
177  if (!h1) {
178  Error("Add","Attempt to add a non-existing profile");
179  return kFALSE;
180  }
181  if (!h1->InheritsFrom(TProfile3D::Class())) {
182  Error("Add","Attempt to add a non-profile2D object");
183  return kFALSE;
184  }
185 
186  return TProfileHelper::Add(this, this, h1, 1, c1);
187 }
188 
189 ////////////////////////////////////////////////////////////////////////////////
190 /// Replace contents of this profile3D by the addition of h1 and h2.
191 ///
192 /// `this = c1*h1 + c2*h2`
193 
195 {
196  if (!h1 || !h2) {
197  Error("Add","Attempt to add a non-existing profile");
198  return kFALSE;
199  }
200  if (!h1->InheritsFrom(TProfile3D::Class())) {
201  Error("Add","Attempt to add a non-profile3D object");
202  return kFALSE;
203  }
204  if (!h2->InheritsFrom(TProfile3D::Class())) {
205  Error("Add","Attempt to add a non-profile3D object");
206  return kFALSE;
207  }
208 
209  return TProfileHelper::Add(this, h1, h2, c1, c2);
210 }
211 
212 
213 ////////////////////////////////////////////////////////////////////////////////
214 /// Set the fgApproximate flag.
215 ///
216 /// When the flag is true, the function GetBinError
217 /// will approximate the bin error with the average profile error on all bins
218 /// in the following situation only
219 ///
220 /// - the number of bins in the profile3D is less than 10404 (eg 100x100x100)
221 /// - the bin number of entries is small ( <5)
222 /// - the estimated bin error is extremely small compared to the bin content
223 /// (see TProfile3D::GetBinError)
224 
226 {
227  fgApproximate = approx;
228 }
229 
230 
231 ////////////////////////////////////////////////////////////////////////////////
232 /// Fill histogram with all entries in the buffer.
233 ///
234 /// - action = -1 histogram is reset and refilled from the buffer (called by THistPainter::Paint)
235 /// - action = 0 histogram is filled from the buffer
236 /// - action = 1 histogram is filled and buffer is deleted
237 /// The buffer is automatically deleted when the number of entries
238 /// in the buffer is greater than the number of entries in the histogram
239 
241 {
242  // do we need to compute the bin size?
243  if (!fBuffer) return 0;
244  Int_t nbentries = (Int_t)fBuffer[0];
245  if (!nbentries) return 0;
246  Double_t *buffer = fBuffer;
247  if (nbentries < 0) {
248  if (action == 0) return 0;
249  nbentries = -nbentries;
250  fBuffer=0;
251  Reset("ICES"); // reset without deleting the functions
252  fBuffer = buffer;
253  }
255  //find min, max of entries in buffer
256  Double_t xmin = fBuffer[2];
257  Double_t xmax = xmin;
258  Double_t ymin = fBuffer[3];
259  Double_t ymax = ymin;
260  Double_t zmin = fBuffer[4];
261  Double_t zmax = zmin;
262  for (Int_t i=1;i<nbentries;i++) {
263  Double_t x = fBuffer[5*i+2];
264  if (x < xmin) xmin = x;
265  if (x > xmax) xmax = x;
266  Double_t y = fBuffer[5*i+3];
267  if (y < ymin) ymin = y;
268  if (y > ymax) ymax = y;
269  Double_t z = fBuffer[5*i+4];
270  if (z < zmin) zmin = z;
271  if (z > zmax) zmax = z;
272  }
273  if (fXaxis.GetXmax() <= fXaxis.GetXmin() || fYaxis.GetXmax() <= fYaxis.GetXmin() || fZaxis.GetXmax() <= fZaxis.GetXmin()) {
274  THLimitsFinder::GetLimitsFinder()->FindGoodLimits(this,xmin,xmax,ymin,ymax,zmin,zmax);
275  } else {
276  fBuffer = 0;
277  Int_t keep = fBufferSize; fBufferSize = 0;
278  if (xmin < fXaxis.GetXmin()) ExtendAxis(xmin,&fXaxis);
279  if (xmax >= fXaxis.GetXmax()) ExtendAxis(xmax,&fXaxis);
280  if (ymin < fYaxis.GetXmin()) ExtendAxis(ymin,&fYaxis);
281  if (ymax >= fYaxis.GetXmax()) ExtendAxis(ymax,&fYaxis);
282  if (zmin < fZaxis.GetXmin()) ExtendAxis(zmin,&fZaxis);
283  if (zmax >= fZaxis.GetXmax()) ExtendAxis(zmax,&fZaxis);
284  fBuffer = buffer;
285  fBufferSize = keep;
286  }
287  }
288 
289  fBuffer = 0;
290  for (Int_t i=0;i<nbentries;i++) {
291  Fill(buffer[5*i+2],buffer[5*i+3],buffer[5*i+4],buffer[5*i+5],buffer[5*i+1]);
292  }
293  fBuffer = buffer;
294 
295  if (action > 0) { delete [] fBuffer; fBuffer = 0; fBufferSize = 0;}
296  else {
297  if (nbentries == (Int_t)fEntries) fBuffer[0] = -nbentries;
298  else fBuffer[0] = 0;
299  }
300  return nbentries;
301 }
302 
303 ////////////////////////////////////////////////////////////////////////////////
304 /// Accumulate arguments in buffer.
305 ///
306 /// When buffer is full, empty the buffer
307 ///
308 /// - fBuffer[0] = number of entries in buffer
309 /// - fBuffer[1] = w of first entry
310 /// - fBuffer[2] = x of first entry
311 /// - fBuffer[3] = y of first entry
312 /// - fBuffer[4] = z of first entry
313 /// - fBuffer[5] = t of first entry
314 
316 {
317  if (!fBuffer) return -3;
318  Int_t nbentries = (Int_t)fBuffer[0];
319  if (nbentries < 0) {
320  nbentries = -nbentries;
321  fBuffer[0] = nbentries;
322  if (fEntries > 0) {
323  Double_t *buffer = fBuffer; fBuffer=0;
324  Reset("ICES"); // reset without deleting the functions
325  fBuffer = buffer;
326  }
327  }
328  if (5*nbentries+5 >= fBufferSize) {
329  BufferEmpty(1);
330  return Fill(x,y,z,t,w);
331  }
332  fBuffer[5*nbentries+1] = w;
333  fBuffer[5*nbentries+2] = x;
334  fBuffer[5*nbentries+3] = y;
335  fBuffer[5*nbentries+4] = z;
336  fBuffer[5*nbentries+5] = t;
337  fBuffer[0] += 1;
338  return -2;
339 }
340 
341 ////////////////////////////////////////////////////////////////////////////////
342 /// Copy a Profile3D histogram to a new profile2D histogram.
343 
344 void TProfile3D::Copy(TObject &obj) const
345 {
346  try {
347  TProfile3D & pobj = dynamic_cast<TProfile3D&>(obj);
348 
349  TH3D::Copy(pobj);
351  fBinSumw2.Copy(pobj.fBinSumw2);
352  for (int bin=0;bin<fNcells;bin++) {
353  pobj.fArray[bin] = fArray[bin];
354  pobj.fSumw2.fArray[bin] = fSumw2.fArray[bin];
355  }
356  pobj.fTmin = fTmin;
357  pobj.fTmax = fTmax;
358  pobj.fScaling = fScaling;
359  pobj.fErrorMode = fErrorMode;
360  pobj.fTsumwt = fTsumwt;
361  pobj.fTsumwt2 = fTsumwt2;
362 
363  } catch(...) {
364  Fatal("Copy","Cannot copy a TProfile3D in a %s",obj.IsA()->GetName());
365  }
366 }
367 
368 ////////////////////////////////////////////////////////////////////////////////
369 /// Performs the operation: `this = this/(c1*f1)` .
370 ///
371 /// This function is not implemented
372 
374 {
375  Error("Divide","Function not implemented for TProfile3D");
376  return kFALSE;
377 }
378 
379 ////////////////////////////////////////////////////////////////////////////////
380 /// Divide this profile2D by h1.
381 ///
382 /// `this = this/h1`
383 ///
384 /// This function return kFALSE if the divide operation failed
385 
387 {
388  if (!h1) {
389  Error("Divide","Attempt to divide a non-existing profile2D");
390  return kFALSE;
391  }
392  if (!h1->InheritsFrom(TProfile3D::Class())) {
393  Error("Divide","Attempt to divide a non-profile3D object");
394  return kFALSE;
395  }
396  TProfile3D *p1 = (TProfile3D*)h1;
397 
398  // delete buffer if it is there since it will become invalid
399  if (fBuffer) BufferEmpty(1);
400 
401 // Check profile compatibility
402  Int_t nx = GetNbinsX();
403  if (nx != p1->GetNbinsX()) {
404  Error("Divide","Attempt to divide profiles with different number of bins");
405  return kFALSE;
406  }
407  Int_t ny = GetNbinsY();
408  if (ny != p1->GetNbinsY()) {
409  Error("Divide","Attempt to divide profiles with different number of bins");
410  return kFALSE;
411  }
412  Int_t nz = GetNbinsZ();
413  if (nz != p1->GetNbinsZ()) {
414  Error("Divide","Attempt to divide profiles with different number of bins");
415  return kFALSE;
416  }
417 
418 // Reset statistics
419  fEntries = fTsumw = fTsumw2 = fTsumwx = fTsumwx2 = 0;
420 
421 // Loop on bins (including underflows/overflows)
422  Int_t bin,binx,biny,binz;
423  Double_t *cu1 = p1->GetW();
424  Double_t *er1 = p1->GetW2();
425  Double_t *en1 = p1->GetB();
426  Double_t c0,c1,w,u,x,y,z;
427  for (binx =0;binx<=nx+1;binx++) {
428  for (biny =0;biny<=ny+1;biny++) {
429  for (binz =0;binz<=nz+1;binz++) {
430  bin = GetBin(binx,biny,binz);
431  c0 = fArray[bin];
432  c1 = cu1[bin];
433  if (c1) w = c0/c1;
434  else w = 0;
435  fArray[bin] = w;
436  u = TMath::Abs(w);
437  x = fXaxis.GetBinCenter(binx);
438  y = fYaxis.GetBinCenter(biny);
439  z = fZaxis.GetBinCenter(binz);
440  fEntries++;
441  fTsumw += u;
442  fTsumw2 += u*u;
443  fTsumwx += u*x;
444  fTsumwx2 += u*x*x;
445  fTsumwy += u*y;
446  fTsumwy2 += u*y*y;
447  fTsumwxy += u*x*y;
448  fTsumwz += u;
449  fTsumwz2 += u*z;
450  fTsumwxz += u*x*z;
451  fTsumwyz += u*y*z;
452  fTsumwt += u;
453  fTsumwt2 += u*u;
454  Double_t e0 = fSumw2.fArray[bin];
455  Double_t e1 = er1[bin];
456  Double_t c12= c1*c1;
457  if (!c1) fSumw2.fArray[bin] = 0;
458  else fSumw2.fArray[bin] = (e0*c1*c1 + e1*c0*c0)/(c12*c12);
459  if (!en1[bin]) fBinEntries.fArray[bin] = 0;
460  else fBinEntries.fArray[bin] /= en1[bin];
461  }
462  }
463  }
464  // maintaining the correct sum of weights square is not supported when dividing
465  // bin error resulting from division of profile needs to be checked
466  if (fBinSumw2.fN) {
467  Warning("Divide","Cannot preserve during the division of profiles the sum of bin weight square");
468  fBinSumw2 = TArrayD();
469  }
470  return kTRUE;
471 }
472 
473 ////////////////////////////////////////////////////////////////////////////////
474 /// Replace contents of this profile2D by the division of h1 by h2.
475 ///
476 /// `this = c1*h1/(c2*h2)`
477 ///
478 /// This function return kFALSE if the divide operation failed
479 
481 {
482  TString opt = option;
483  opt.ToLower();
484  Bool_t binomial = kFALSE;
485  if (opt.Contains("b")) binomial = kTRUE;
486  if (!h1 || !h2) {
487  Error("Divide","Attempt to divide a non-existing profile2D");
488  return kFALSE;
489  }
490  if (!h1->InheritsFrom(TProfile3D::Class())) {
491  Error("Divide","Attempt to divide a non-profile2D object");
492  return kFALSE;
493  }
494  TProfile3D *p1 = (TProfile3D*)h1;
495  if (!h2->InheritsFrom(TProfile3D::Class())) {
496  Error("Divide","Attempt to divide a non-profile2D object");
497  return kFALSE;
498  }
499  TProfile3D *p2 = (TProfile3D*)h2;
500 
501 // Check histogram compatibility
502  Int_t nx = GetNbinsX();
503  if (nx != p1->GetNbinsX() || nx != p2->GetNbinsX()) {
504  Error("Divide","Attempt to divide profiles with different number of bins");
505  return kFALSE;
506  }
507  Int_t ny = GetNbinsY();
508  if (ny != p1->GetNbinsY() || ny != p2->GetNbinsY()) {
509  Error("Divide","Attempt to divide profiles with different number of bins");
510  return kFALSE;
511  }
512  Int_t nz = GetNbinsZ();
513  if (nz != p1->GetNbinsZ() || nz != p2->GetNbinsZ()) {
514  Error("Divide","Attempt to divide profiles with different number of bins");
515  return kFALSE;
516  }
517  if (!c2) {
518  Error("Divide","Coefficient of dividing profile cannot be zero");
519  return kFALSE;
520  }
521 
522 // Reset statistics
523  fEntries = fTsumw = fTsumw2 = fTsumwx = fTsumwx2 = 0;
524 
525 // Loop on bins (including underflows/overflows)
526  Int_t bin,binx,biny,binz;
527  Double_t *cu1 = p1->GetW();
528  Double_t *cu2 = p2->GetW();
529  Double_t *er1 = p1->GetW2();
530  Double_t *er2 = p2->GetW2();
531  Double_t *en1 = p1->GetB();
532  Double_t *en2 = p2->GetB();
533  Double_t b1,b2,w,u,x,y,z,ac1,ac2;
534  ac1 = TMath::Abs(c1);
535  ac2 = TMath::Abs(c2);
536  for (binx =0;binx<=nx+1;binx++) {
537  for (biny =0;biny<=ny+1;biny++) {
538  for (binz =0;binz<=nz+1;binz++) {
539  bin = GetBin(binx,biny,binz);
540  b1 = cu1[bin];
541  b2 = cu2[bin];
542  if (b2) w = c1*b1/(c2*b2);
543  else w = 0;
544  fArray[bin] = w;
545  u = TMath::Abs(w);
546  x = fXaxis.GetBinCenter(binx);
547  y = fYaxis.GetBinCenter(biny);
548  z = fZaxis.GetBinCenter(biny);
549  fEntries++;
550  fTsumw += u;
551  fTsumw2 += u*u;
552  fTsumwx += u*x;
553  fTsumwx2 += u*x*x;
554  fTsumwy += u*y;
555  fTsumwy2 += u*y*y;
556  fTsumwxy += u*x*y;
557  fTsumwz += u*z;
558  fTsumwz2 += u*z*z;
559  fTsumwxz += u*x*z;
560  fTsumwyz += u*y*z;
561  fTsumwt += u;
562  fTsumwt2 += u*u;
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  }
580  return kTRUE;
581 }
582 
583 ////////////////////////////////////////////////////////////////////////////////
584 /// Fill a Profile3D histogram (no weights).
585 
587 {
588  if (fBuffer) return BufferFill(x,y,z,t,1);
589 
590  Int_t bin,binx,biny,binz;
591 
592  if (fTmin != fTmax) {
593  if (t <fTmin || t> fTmax || TMath::IsNaN(t) ) return -1;
594  }
595 
596  fEntries++;
597  binx =fXaxis.FindBin(x);
598  biny =fYaxis.FindBin(y);
599  binz =fZaxis.FindBin(z);
600  if (binx <0 || biny <0 || binz<0) return -1;
601  bin = GetBin(binx,biny,binz);
602  AddBinContent(bin, t);
603  fSumw2.fArray[bin] += (Double_t)t*t;
604  fBinEntries.fArray[bin] += 1;
605  if (fBinSumw2.fN) fBinSumw2.fArray[bin] += 1;
606  if (binx == 0 || binx > fXaxis.GetNbins()) {
607  if (!fgStatOverflows) return -1;
608  }
609  if (biny == 0 || biny > fYaxis.GetNbins()) {
610  if (!fgStatOverflows) return -1;
611  }
612  if (binz == 0 || binz > fZaxis.GetNbins()) {
613  if (!fgStatOverflows) return -1;
614  }
615 //printf("x=%g, y=%g, z=%g, t=%g, binx=%d, biny=%d, binz=%d, bin=%d\n",x,y,z,t,binx,biny,binz,bin);
616  ++fTsumw;
617  ++fTsumw2;
618  fTsumwx += x;
619  fTsumwx2 += x*x;
620  fTsumwy += y;
621  fTsumwy2 += y*y;
622  fTsumwxy += x*y;
623  fTsumwz += z;
624  fTsumwz2 += z*z;
625  fTsumwxz += x*z;
626  fTsumwyz += y*z;
627  fTsumwt += t;
628  fTsumwt2 += t*t;
629  return bin;
630 }
631 
632 ////////////////////////////////////////////////////////////////////////////////
633 /// Fill a Profile3D histogram with weights.
634 
636 {
637  if (fBuffer) return BufferFill(x,y,z,t,w);
638 
639  Int_t bin,binx,biny,binz;
640 
641  if (fTmin != fTmax) {
642  if (t <fTmin || t> fTmax || TMath::IsNaN(t) ) return -1;
643  }
644 
645  Double_t u= w; // (w > 0 ? w : -w);
646  fEntries++;
647  binx =fXaxis.FindBin(x);
648  biny =fYaxis.FindBin(y);
649  binz =fZaxis.FindBin(z);
650  if (binx <0 || biny <0 || binz<0) return -1;
651  bin = GetBin(binx,biny,binz);
652  AddBinContent(bin, u*t);
653  fSumw2.fArray[bin] += u*t*t;
654  if (!fBinSumw2.fN && u != 1.0 && !TestBit(TH1::kIsNotW)) Sumw2(); // must be called before accumulating the entries
655  if (fBinSumw2.fN) fBinSumw2.fArray[bin] += u*u;
656  fBinEntries.fArray[bin] += u;
657  if (binx == 0 || binx > fXaxis.GetNbins()) {
658  if (!fgStatOverflows) return -1;
659  }
660  if (biny == 0 || biny > fYaxis.GetNbins()) {
661  if (!fgStatOverflows) return -1;
662  }
663  if (binz == 0 || binz > fZaxis.GetNbins()) {
664  if (!fgStatOverflows) return -1;
665  }
666  fTsumw += u;
667  fTsumw2 += u*u;
668  fTsumwx += u*x;
669  fTsumwx2 += u*x*x;
670  fTsumwy += u*y;
671  fTsumwy2 += u*y*y;
672  fTsumwxy += u*x*y;
673  fTsumwz += u*z;
674  fTsumwz2 += u*z*z;
675  fTsumwxz += u*x*z;
676  fTsumwyz += u*y*z;
677  fTsumwt += u*t;
678  fTsumwt2 += u*t*t;
679  return bin;
680 }
681 
682 ////////////////////////////////////////////////////////////////////////////////
683 /// Return bin content of a Profile3D histogram.
684 
686 {
687  if (fBuffer) ((TProfile3D*)this)->BufferEmpty();
688 
689  if (bin < 0 || bin >= fNcells) return 0;
690  if (fBinEntries.fArray[bin] == 0) return 0;
691  if (!fArray) return 0;
692  return fArray[bin]/fBinEntries.fArray[bin];
693 }
694 
695 ////////////////////////////////////////////////////////////////////////////////
696 /// Return bin entries of a Profile3D histogram.
697 
699 {
700  if (fBuffer) ((TProfile3D*)this)->BufferEmpty();
701 
702  if (bin < 0 || bin >= fNcells) return 0;
703  return fBinEntries.fArray[bin];
704 }
705 
706 ////////////////////////////////////////////////////////////////////////////////
707 /// Return bin effective entries for a weighted filled Profile histogram.
708 ///
709 /// In case of an unweighted profile, it is equivalent to the number of entries per bin
710 /// The effective entries is defined as the square of the sum of the weights divided by the
711 /// sum of the weights square.
712 /// TProfile::Sumw2() must be called before filling the profile with weights.
713 /// Only by calling this method the sum of the square of the weights per bin is stored.
714 
716 {
718 }
719 
720 ////////////////////////////////////////////////////////////////////////////////
721 /// Return bin error of a Profile3D histogram.
722 ///
723 /// ### Computing errors: A moving field
724 ///
725 /// The computation of errors for a TProfile3D has evolved with the versions
726 /// of ROOT. The difficulty is in computing errors for bins with low statistics.
727 ///
728 /// - prior to version 3.10, we had no special treatment of low statistic bins.
729 /// As a result, these bins had huge errors. The reason is that the
730 /// expression eprim2 is very close to 0 (rounding problems) or 0.
731 /// - The algorithm is modified/protected for the case
732 /// when a TProfile3D is projected (ProjectionX). The previous algorithm
733 /// generated a N^2 problem when projecting a TProfile3D with a large number of
734 /// bins (eg 100000).
735 /// - in version 3.10/02, a new static function TProfile::Approximate
736 /// is introduced to enable or disable (default) the approximation.
737 /// (see also comments in TProfile::GetBinError)
738 
740 {
741  return TProfileHelper::GetBinError((TProfile3D*)this, bin);
742 }
743 
744 ////////////////////////////////////////////////////////////////////////////////
745 /// Return option to compute profile2D errors.
746 
748 {
749  if (fErrorMode == kERRORSPREAD) return "s";
750  if (fErrorMode == kERRORSPREADI) return "i";
751  if (fErrorMode == kERRORSPREADG) return "g";
752  return "";
753 }
754 
755 ////////////////////////////////////////////////////////////////////////////////
756 /// fill the array stats from the contents of this profile.
757 ///
758 /// The array stats must be correctly dimensionned in the calling program.
759 ///
760 /// - stats[0] = sumw
761 /// - stats[1] = sumw2
762 /// - stats[2] = sumwx
763 /// - stats[3] = sumwx2
764 /// - stats[4] = sumwy
765 /// - stats[5] = sumwy2
766 /// - stats[6] = sumwxy
767 /// - stats[7] = sumwz
768 /// - stats[8] = sumwz2
769 /// - stats[9] = sumwxz
770 /// - stats[10]= sumwyz
771 /// - stats[11]= sumwt
772 /// - stats[12]= sumwt2
773 ///
774 /// If no axis-subrange is specified (via TAxis::SetRange), the array stats
775 /// is simply a copy of the statistics quantities computed at filling time.
776 /// If a sub-range is specified, the function recomputes these quantities
777 /// from the bin contents in the current axis range.
778 
779 void TProfile3D::GetStats(Double_t *stats) const
780 {
781  if (fBuffer) ((TProfile3D*)this)->BufferEmpty();
782 
783  // Loop on bins
785  Int_t bin, binx, biny,binz;
786  Double_t w, w2;
787  Double_t x,y,z;
788  for (bin=0;bin<kNstat;bin++) stats[bin] = 0;
789  if (!fBinEntries.fArray) return;
790  for (binz=fZaxis.GetFirst();binz<=fZaxis.GetLast();binz++) {
791  z = fZaxis.GetBinCenter(binz);
792  for (biny=fYaxis.GetFirst();biny<=fYaxis.GetLast();biny++) {
793  y = fYaxis.GetBinCenter(biny);
794  for (binx=fXaxis.GetFirst();binx<=fXaxis.GetLast();binx++) {
795  bin = GetBin(binx,biny,binz);
796  w = fBinEntries.fArray[bin];
797  w2 = (fBinSumw2.fN ? fBinSumw2.fArray[bin] : w );
798  x = fXaxis.GetBinCenter(binx);
799  stats[0] += w;
800  stats[1] += w2;
801  stats[2] += w*x;
802  stats[3] += w*x*x;
803  stats[4] += w*y;
804  stats[5] += w*y*y;
805  stats[6] += w*x*y;
806  stats[7] += w*z;
807  stats[8] += w*z*z;
808  stats[9] += w*x*z;
809  stats[10] += w*y*z;
810  stats[11] += fArray[bin];
811  stats[12] += fSumw2.fArray[bin];
812  }
813  }
814  }
815  } else {
816  stats[0] = fTsumw;
817  stats[1] = fTsumw2;
818  stats[2] = fTsumwx;
819  stats[3] = fTsumwx2;
820  stats[4] = fTsumwy;
821  stats[5] = fTsumwy2;
822  stats[6] = fTsumwxy;
823  stats[7] = fTsumwz;
824  stats[8] = fTsumwz2;
825  stats[9] = fTsumwxz;
826  stats[10] = fTsumwyz;
827  stats[11] = fTsumwt;
828  stats[12] = fTsumwt2;
829  }
830 }
831 
832 ////////////////////////////////////////////////////////////////////////////////
833 /// Merge all histograms in the collection in this histogram.
834 ///
835 /// This function computes the min/max for the axes,
836 /// compute a new number of bins, if necessary,
837 /// add bin contents, errors and statistics.
838 /// If overflows are present and limits are different the function will fail.
839 /// The function returns the total number of entries in the result histogram
840 /// if the merge is successful, -1 otherwise.
841 ///
842 /// IMPORTANT remark. The 2 axis x and y may have different number
843 /// of bins and different limits, BUT the largest bin width must be
844 /// a multiple of the smallest bin width and the upper limit must also
845 /// be a multiple of the bin width.
846 
848 {
849  return TProfileHelper::Merge(this, li);
850 }
851 
852 ////////////////////////////////////////////////////////////////////////////////
853 /// Performs the operation: `this = this*c1*f1` .
854 
856 {
857  Error("Multiply","Function not implemented for TProfile3D");
858  return kFALSE;
859 }
860 
861 ////////////////////////////////////////////////////////////////////////////////
862 /// Multiply this profile2D by h1.
863 ///
864 /// `this = this*h1`
865 
867 {
868  Error("Multiply","Multiplication of profile2D histograms not implemented");
869  return kFALSE;
870 }
871 
872 ////////////////////////////////////////////////////////////////////////////////
873 /// Replace contents of this profile2D by multiplication of h1 by h2.
874 ///
875 /// `this = (c1*h1)*(c2*h2)`
876 
878 {
879  Error("Multiply","Multiplication of profile2D histograms not implemented");
880  return kFALSE;
881 }
882 
883 ////////////////////////////////////////////////////////////////////////////////
884 /// Project this profile3D into a 3-D histogram along X,Y,Z.
885 ///
886 /// The projection is always of the type TH3D.
887 ///
888 /// - if option "E" is specified, the errors are computed. (default)
889 /// - if option "B" is specified, the content of bin of the returned histogram
890 /// will be equal to the GetBinEntries(bin) of the profile,
891 /// - if option "C=E" the bin contents of the projection are set to the
892 /// bin errors of the profile
893 /// - if option "E" is specified the errors of the projected histogram are computed and set
894 /// to be equal to the errors of the profile.
895 /// Option "E" is defined as the default one in the header file.
896 /// - if option "" is specified the histogram errors are simply the sqrt of its content
897 /// - if option "B" is specified, the content of bin of the returned histogram
898 /// will be equal to the GetBinEntries(bin) of the profile,
899 /// - if option "C=E" the bin contents of the projection are set to the
900 /// bin errors of the profile
901 /// - if option "W" is specified the bin content of the projected histogram is set to the
902 /// product of the bin content of the profile and the entries.
903 /// With this option the returned histogram will be equivalent to the one obtained by
904 /// filling directly a TH2D using the 3-rd value as a weight.
905 /// This option makes sense only for profile filled with all weights =1.
906 /// When the profile is weighted (filled with weights different than 1) the
907 /// bin error of the projected histogram (obtained using this option "W") cannot be
908 /// correctly computed from the information stored in the profile. In that case the
909 /// obtained histogram contains as bin error square the weighted sum of the square of the
910 /// profiled observable (TProfile2D::fSumw2[bin] )
911 
912 TH3D *TProfile3D::ProjectionXYZ(const char *name, Option_t *option) const
913 {
914 
915  TString opt = option;
916  opt.ToLower();
917  Int_t nx = fXaxis.GetNbins();
918  Int_t ny = fYaxis.GetNbins();
919  Int_t nz = fZaxis.GetNbins();
920  const TArrayD *xbins = fXaxis.GetXbins();
921  const TArrayD *ybins = fYaxis.GetXbins();
922  const TArrayD *zbins = fZaxis.GetXbins();
923 
924  // Create the projection histogram
925  TString pname = name;
926  if (pname == "_px") {
927  pname = GetName(); pname.Append("_pxyz");
928  }
929  TH3D *h1 = 0 ;
930  if (xbins->fN == 0 && ybins->fN == 0 && zbins->fN == 0)
932  else if ( xbins->fN != 0 && ybins->fN != 0 && zbins->fN != 0)
933  h1 = new TH3D(pname,GetTitle(),nx,xbins->GetArray(),ny,ybins->GetArray(), nz,zbins->GetArray() );
934  else {
935  Error("ProjectionXYZ","Histogram has an axis with variable bins and an axis with fixed bins. This case is not supported - return a null pointer");
936  return 0;
937  }
938 
939 
940  Bool_t computeErrors = kFALSE;
941  Bool_t cequalErrors = kFALSE;
942  Bool_t binEntries = kFALSE;
943  Bool_t binWeight = kFALSE;
944 
945  if (opt.Contains("b")) binEntries = kTRUE;
946  if (opt.Contains("e")) computeErrors = kTRUE;
947  if (opt.Contains("w")) binWeight = kTRUE;
948  if (opt.Contains("c=e")) {cequalErrors = kTRUE; computeErrors=kFALSE;}
949  if (computeErrors || binWeight || (binEntries && fBinSumw2.fN) ) h1->Sumw2();
950 
951  // Fill the projected histogram
952  Int_t bin,binx,biny,binz;
953  Double_t cont;
954  for (binx =0;binx<=nx+1;binx++) {
955  for (biny =0;biny<=ny+1;biny++) {
956  for (binz =0;binz<=nz+1;binz++) {
957  bin = GetBin(binx,biny,binz);
958 
959  if (binEntries) cont = GetBinEntries(bin);
960  else if (cequalErrors) cont = GetBinError(bin);
961  else if (binWeight) cont = GetBinContent(bin) * GetBinEntries(bin);
962  else cont = GetBinContent(bin); // default case
963 
964  h1->SetBinContent(bin ,cont);
965 
966  // if option E projected histogram errors are same as profile
967  if (computeErrors ) h1->SetBinError(bin , GetBinError(bin) );
968  // in case of option W bin error is deduced from bin sum of z**2 values of profile
969  // this is correct only if the profile is unweighted
970  if (binWeight) h1->GetSumw2()->fArray[bin] = fSumw2.fArray[bin];
971  // in case of bin entries and profile is weighted, we need to set also the bin error
972  if (binEntries && fBinSumw2.fN ) {
973  R__ASSERT( h1->GetSumw2() );
974  h1->GetSumw2()->fArray[bin] = fBinSumw2.fArray[bin];
975  }
976  }
977  }
978  }
979  h1->SetEntries(fEntries);
980  return h1;
981 }
982 ////////////////////////////////////////////////////////////////////////////////
983 /// Project a 3-D profile into a 2D-profile histogram depending on the option parameter.
984 ///
985 /// option may contain a combination of the characters x,y,z:
986 ///
987 /// - option = "xy" return the x versus y projection into a TProfile2D histogram
988 /// - option = "yx" return the y versus x projection into a TProfile2D histogram
989 /// - option = "xz" return the x versus z projection into a TProfile2D histogram
990 /// - option = "zx" return the z versus x projection into a TProfile2D histogram
991 /// - option = "yz" return the y versus z projection into a TProfile2D histogram
992 /// - option = "zy" return the z versus y projection into a TProfile2D histogram
993 ///
994 /// NB: the notation "a vs b" means "a" vertical and "b" horizontal along X
995 ///
996 /// The resulting profile contains the combination of all the considered bins along X
997 /// By default, all bins are included considering also underflow/overflows
998 ///
999 /// The option can also be used to specify the projected profile error type.
1000 /// Values which can be used are 's', 'i', or 'g'. See TProfile::BuildOptions for details
1001 ///
1002 /// To select a bin range along an axis, use TAxis::SetRange, eg
1003 /// `h3.GetYaxis()->SetRange(23,56);`
1004 
1006 {
1007  // can call TH3 method which will call the virtual method :DoProjectProfile2D re-implemented below
1008  // but need to add underflow/overflow
1009  TString opt(option);
1010  opt.Append(" UF OF");
1011  return TH3::Project3DProfile(opt);
1012 }
1013 
1014 ////////////////////////////////////////////////////////////////////////////////
1015 /// Internal method to project to a 2D Profile.
1016 ///
1017 /// Called from TH3::Project3DProfile but re-implemented in case of the TPRofile3D
1018 /// since what is done is different.
1019 
1020 TProfile2D *TProfile3D::DoProjectProfile2D(const char* name, const char * title, const TAxis* projX, const TAxis* projY,
1021  bool originalRange, bool useUF, bool useOF) const
1022 {
1023  // Get the ranges where we will work.
1024  Int_t ixmin = projX->GetFirst();
1025  Int_t ixmax = projX->GetLast();
1026  Int_t iymin = projY->GetFirst();
1027  Int_t iymax = projY->GetLast();
1028  if (ixmin == 0 && ixmax == 0) { ixmin = 1; ixmax = projX->GetNbins(); }
1029  if (iymin == 0 && iymax == 0) { iymin = 1; iymax = projY->GetNbins(); }
1030  Int_t nx = ixmax-ixmin+1;
1031  Int_t ny = iymax-iymin+1;
1032 
1033  // Create the projected profiles
1034  TProfile2D *p2 = 0;
1035  // Create always a new TProfile2D (not as in the case of TH3 projection)
1036 
1037  const TArrayD *xbins = projX->GetXbins();
1038  const TArrayD *ybins = projY->GetXbins();
1039  // assume all axis have variable bins or have fixed bins
1040  if ( originalRange ) {
1041  if (xbins->fN == 0 && ybins->fN == 0) {
1042  p2 = new TProfile2D(name,title,projY->GetNbins(),projY->GetXmin(),projY->GetXmax()
1043  ,projX->GetNbins(),projX->GetXmin(),projX->GetXmax());
1044  } else {
1045  p2 = new TProfile2D(name,title,projY->GetNbins(),&ybins->fArray[iymin-1],projX->GetNbins(),&xbins->fArray[ixmin-1]);
1046  }
1047  } else {
1048  if (xbins->fN == 0 && ybins->fN == 0) {
1049  p2 = new TProfile2D(name,title,ny,projY->GetBinLowEdge(iymin),projY->GetBinUpEdge(iymax)
1050  ,nx,projX->GetBinLowEdge(ixmin),projX->GetBinUpEdge(ixmax));
1051  } else {
1052  p2 = new TProfile2D(name,title,ny,&ybins->fArray[iymin-1],nx,&xbins->fArray[ixmin-1]);
1053  }
1054  }
1055 
1056  // weights
1057  bool useWeights = (fBinSumw2.fN != 0);
1058  if (useWeights) p2->Sumw2();
1059 
1060  // make projection in a 3D first
1061  TH3D * h3dW = ProjectionXYZ("h3temp-W","W");
1062  TH3D * h3dN = ProjectionXYZ("h3temp-N","B");
1063 
1064  h3dW->SetDirectory(0); h3dN->SetDirectory(0);
1065 
1066  // note that h3dW is always a weighted histogram - so we need to compute error in the projection
1067  TAxis * projX_hW = h3dW->GetXaxis();
1068  TAxis * projX_hN = h3dN->GetXaxis();
1069  if (projX == GetYaxis() ) { projX_hW = h3dW->GetYaxis(); projX_hN = h3dN->GetYaxis(); }
1070  if (projX == GetZaxis() ) { projX_hW = h3dW->GetZaxis(); projX_hN = h3dN->GetZaxis(); }
1071  TAxis * projY_hW = h3dW->GetYaxis();
1072  TAxis * projY_hN = h3dN->GetYaxis();
1073  if (projY == GetXaxis() ) { projY_hW = h3dW->GetXaxis(); projY_hN = h3dN->GetXaxis(); }
1074  if (projY == GetZaxis() ) { projY_hW = h3dW->GetZaxis(); projY_hN = h3dN->GetZaxis(); }
1075 
1076  TH2D * h2W = TH3::DoProject2D(*h3dW,"htemp-W","",projX_hW, projY_hW, true, originalRange, useUF, useOF);
1077  TH2D * h2N = TH3::DoProject2D(*h3dN,"htemp-N","",projX_hN, projY_hN, useWeights, originalRange, useUF, useOF);
1078  h2W->SetDirectory(0); h2N->SetDirectory(0);
1079 
1080 
1081  // fill the bin content
1082  R__ASSERT( h2W->fN == p2->fN );
1083  R__ASSERT( h2N->fN == p2->fN );
1084  R__ASSERT( h2W->GetSumw2()->fN != 0); // h2W should always be a weighted histogram since h3dW is weighted
1085  for (int i = 0; i < p2->fN ; ++i) {
1086  //std::cout << " proj bin " << i << " " << h2W->fArray[i] << " " << h2N->fArray[i] << std::endl;
1087  p2->fArray[i] = h2W->fArray[i]; // array of profile is sum of all values
1088  p2->GetSumw2()->fArray[i] = h2W->GetSumw2()->fArray[i]; // array of content square of profile is weight square of the W projected histogram
1089  p2->SetBinEntries(i, h2N->fArray[i] );
1090  if (useWeights) p2->GetBinSumw2()->fArray[i] = h2N->GetSumw2()->fArray[i]; // sum of weight squares are stored to compute errors in h1N histogram
1091  }
1092  // delete the created histograms
1093  delete h3dW;
1094  delete h3dN;
1095  delete h2W;
1096  delete h2N;
1097 
1098  // Also we need to set the entries since they have not been correctly calculated during the projection
1099  // we can only set them to the effective entries
1100  p2->SetEntries( p2->GetEffectiveEntries() );
1101 
1102  return p2;
1103 
1104 }
1105 
1106 ////////////////////////////////////////////////////////////////////////////////
1107 /// Replace current statistics with the values in array stats.
1108 
1110 {
1111  TH3::PutStats(stats);
1112  fTsumwt = stats[11];
1113  fTsumwt2 = stats[12];
1114 }
1115 
1116 ////////////////////////////////////////////////////////////////////////////////
1117 /// Reset contents of a Profile3D histogram.
1118 
1120 {
1121  TH3D::Reset(option);
1122  fBinSumw2.Reset();
1123  fBinEntries.Reset();
1124  TString opt = option;
1125  opt.ToUpper();
1126  if (opt.Contains("ICE") && !opt.Contains("S")) return;
1127  fTsumwt = fTsumwt2 = 0;
1128 }
1129 
1130 ////////////////////////////////////////////////////////////////////////////////
1131 /// Profile histogram is resized along axis such that x is in the axis range.
1132 /// The new axis limits are recomputed by doubling iteratively
1133 /// the current axis range until the specified value x is within the limits.
1134 /// The algorithm makes a copy of the histogram, then loops on all bins
1135 /// of the old histogram to fill the rebinned histogram.
1136 /// Takes into account errors (Sumw2) if any.
1137 /// The axis must be rebinnable before invoking this function.
1138 /// Ex: `h->GetXaxis()->SetCanExtend(kTRUE)`
1139 
1141 {
1142  TProfile3D* hold = TProfileHelper::ExtendAxis(this, x, axis);
1143  if ( hold ) {
1144  fTsumwt = hold->fTsumwt;
1145  fTsumwt2 = hold->fTsumwt2;
1146  delete hold;
1147  }
1148 }
1149 
1150 ////////////////////////////////////////////////////////////////////////////////
1151 /// Save primitive as a C++ statement(s) on output stream out.
1152 ///
1153 /// Note the following restrictions in the code generated:
1154 /// - variable bin size not implemented
1155 /// - SetErrorOption not implemented
1156 
1157 void TProfile3D::SavePrimitive(std::ostream &out, Option_t *option /*= ""*/)
1158 {
1159  char quote = '"';
1160  out <<" "<<std::endl;
1161  out <<" "<<ClassName()<<" *";
1162 
1163  out << GetName() << " = new " << ClassName() << "(" << quote
1164  << GetName() << quote << "," << quote<< GetTitle() << quote
1165  << "," << GetXaxis()->GetNbins();
1166  out << "," << GetXaxis()->GetXmin()
1167  << "," << GetXaxis()->GetXmax();
1168  out << "," << GetYaxis()->GetNbins();
1169  out << "," << GetYaxis()->GetXmin()
1170  << "," << GetYaxis()->GetXmax();
1171  out << "," << GetZaxis()->GetNbins();
1172  out << "," << GetZaxis()->GetXmin()
1173  << "," << GetZaxis()->GetXmax();
1174  out << "," << fTmin
1175  << "," << fTmax;
1176  out << ");" << std::endl;
1177 
1178 
1179  // save bin entries
1180  Int_t bin;
1181  for (bin=0;bin<fNcells;bin++) {
1182  Double_t bi = GetBinEntries(bin);
1183  if (bi) {
1184  out<<" "<<GetName()<<"->SetBinEntries("<<bin<<","<<bi<<");"<<std::endl;
1185  }
1186  }
1187  //save bin contents
1188  for (bin=0;bin<fNcells;bin++) {
1189  Double_t bc = fArray[bin];
1190  if (bc) {
1191  out<<" "<<GetName()<<"->SetBinContent("<<bin<<","<<bc<<");"<<std::endl;
1192  }
1193  }
1194  // save bin errors
1195  if (fSumw2.fN) {
1196  for (bin=0;bin<fNcells;bin++) {
1197  Double_t be = TMath::Sqrt(fSumw2.fArray[bin]);
1198  if (be) {
1199  out<<" "<<GetName()<<"->SetBinError("<<bin<<","<<be<<");"<<std::endl;
1200  }
1201  }
1202  }
1203 
1204  TH1::SavePrimitiveHelp(out, GetName(), option);
1205 }
1206 
1207 ////////////////////////////////////////////////////////////////////////////////
1208 /// Multiply this profile2D by a constant c1.
1209 ///
1210 /// `this = c1*this`
1211 ///
1212 /// This function uses the services of TProfile3D::Add
1213 
1215 {
1216  TProfileHelper::Scale(this, c1, option);
1217 }
1218 
1219 ////////////////////////////////////////////////////////////////////////////////
1220 ///Set the number of entries in bin.
1221 
1223 {
1224  TProfileHelper::SetBinEntries(this, bin, w);
1225 }
1226 
1227 ////////////////////////////////////////////////////////////////////////////////
1228 /// Redefine x, y and z axis parameters.
1229 
1231 {
1232  TH1::SetBins(nx, xmin, xmax, ny, ymin, ymax, nz, zmin, zmax);
1235 }
1236 
1237 ////////////////////////////////////////////////////////////////////////////////
1238 /// Redefine x, y and z axis parameters with variable bin sizes
1239 
1240 void TProfile3D::SetBins(Int_t nx, const Double_t *xBins, Int_t ny, const Double_t *yBins, Int_t nz, const Double_t *zBins)
1241 {
1242  TH1::SetBins(nx,xBins,ny,yBins,nz,zBins);
1245 }
1246 
1247 ////////////////////////////////////////////////////////////////////////////////
1248 /// Set total number of bins including under/overflow.
1249 ///
1250 /// Reallocate bin contents array
1251 
1253 {
1256 }
1257 
1258 ////////////////////////////////////////////////////////////////////////////////
1259 /// Set the buffer size in units of 8 bytes (double).
1260 
1262 {
1263  if (fBuffer) {
1264  BufferEmpty();
1265  delete [] fBuffer;
1266  fBuffer = 0;
1267  }
1268  if (buffersize <= 0) {
1269  fBufferSize = 0;
1270  return;
1271  }
1272  if (buffersize < 100) buffersize = 100;
1273  fBufferSize = 1 + 5*buffersize;
1274  fBuffer = new Double_t[fBufferSize];
1275  memset(fBuffer,0,sizeof(Double_t)*fBufferSize);
1276 }
1277 
1278 ////////////////////////////////////////////////////////////////////////////////
1279 /// Set option to compute profile3D errors.
1280 ///
1281 /// The computation of the bin errors is based on the parameter option:
1282 /// - ' ' (Default) The bin errors are the standard error on the mean of the bin profiled values (T),
1283 /// i.e. the standard error of the bin contents.
1284 /// Note that if TProfile3D::Approximate() is called, an approximation is used when
1285 /// the spread in T is 0 and the number of bin entries is > 0
1286 /// - 's' The bin errors are the standard deviations of the T bin values
1287 /// Note that if TProfile3D::Approximate() is called, an approximation is used when
1288 /// the spread in T is 0 and the number of bin entries is > 0
1289 /// - 'i' Errors are as in default case (standard errors of the bin contents)
1290 /// The only difference is for the case when the spread in T is zero.
1291 /// In this case for N > 0 the error is 1./SQRT(12.*N)
1292 /// - 'g' Errors are 1./SQRT(W) for W not equal to 0 and 0 for W = 0.
1293 /// W is the sum in the bin of the weights of the profile.
1294 /// This option is for combining measurements t +/- dt,
1295 /// and the profile is filled with values t and weights w = 1/dt**2
1296 ///
1297 /// See TProfile::BuildOptions for explanation of all options
1298 
1300 {
1301  TProfileHelper::SetErrorOption(this, option);
1302 }
1303 
1304 ////////////////////////////////////////////////////////////////////////////////
1305 /// Create/Delete structure to store sum of squares of weights per bin
1306 /// This is needed to compute the correct statistical quantities
1307 /// of a profile filled with weights
1308 ///
1309 /// This function is automatically called when the histogram is created
1310 /// if the static function TH1::SetDefaultSumw2 has been called before.
1311 /// If flag = false the structure is deleted
1312 
1314 {
1315  TProfileHelper::Sumw2(this, flag);
1316 }
const int nx
Definition: kalman.C:16
TArrayD()
Default TArrayD ctor.
Definition: TArrayD.cxx:26
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
static void SetBinEntries(T *p, Int_t bin, Double_t w)
Double_t fTmin
Definition: TProfile3D.h:35
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
Double_t fTsumwt2
Definition: TProfile3D.h:39
virtual Double_t GetBinError(Int_t bin) const
Return bin error of a Profile3D histogram.
Definition: TProfile3D.cxx:739
virtual void SetBinEntries(Int_t bin, Double_t w)
Set the number of entries in bin.
virtual void SetBuffer(Int_t buffersize, Option_t *opt="")
Set the buffer size in units of 8 bytes (double).
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 TH3D * ProjectionXYZ(const char *name="_pxyz", Option_t *option="e") const
Project this profile3D into a 3-D histogram along X,Y,Z.
Definition: TProfile3D.cxx:912
virtual TProfile2D * DoProjectProfile2D(const char *name, const char *title, const TAxis *projX, const TAxis *projY, bool originalRange, bool useUF, bool useOF) const
Internal method to project to a 2D Profile.
virtual void AddBinContent(Int_t bin)
Increment bin content by 1.
Definition: TH3.h:318
virtual void Sumw2(Bool_t flag=kTRUE)
Create/Delete structure to store sum of squares of weights per bin.
virtual void Scale(Double_t c1=1, Option_t *option="")
Multiply this profile2D by a constant c1.
static Double_t GetBinError(T *p, Int_t bin)
Option_t * GetErrorOption() const
Return option to compute profile2D errors.
Definition: TProfile3D.cxx:747
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
virtual void Sumw2(Bool_t flag=kTRUE)
Create/Delete structure to store sum of squares of weights per bin This is needed to compute the corr...
void Reset()
Definition: TArrayD.h:47
float ymin
Definition: THbookFile.cxx:93
TAxis fYaxis
Y axis descriptor.
Definition: TH1.h:81
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content.
Definition: TH3.cxx:3226
static void SetErrorOption(T *p, Option_t *opt)
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
virtual ~TProfile3D()
Default destructor for Profile3D histograms.
Definition: TProfile3D.cxx:87
Double_t fTmax
Definition: TProfile3D.h:36
static Bool_t fgStatOverflows
!flag to use under/overflows in statistics
Definition: TH1.h:106
virtual Double_t GetBinEntries(Int_t bin) const
Return bin entries of a Profile3D histogram.
Definition: TProfile3D.cxx:698
virtual Int_t BufferFill(Double_t, Double_t)
accumulate arguments in buffer.
Definition: TProfile3D.h:43
virtual Int_t GetBin(Int_t binx, Int_t biny, Int_t binz) const
See comments in TH1::GetBin.
Definition: TH3.cxx:950
void ToUpper()
Change string to upper case.
Definition: TString.cxx:1112
virtual Int_t GetNbinsZ() const
Definition: TH1.h:279
#define R__ASSERT(e)
Definition: TError.h:96
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save primitive as a C++ statement(s) on output stream out.
static THLimitsFinder * GetLimitsFinder()
Return pointer to the current finder.
Double_t fTsumwy
Definition: TH3.h:34
Basic string class.
Definition: TString.h:129
void BuildOptions(Double_t tmin, Double_t tmax, Option_t *option)
Set Profile3D histogram structure and options.
Definition: TProfile3D.cxx:142
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
Profile3D histograms are used to display the mean value of T and its RMS for each cell in X...
Definition: TProfile3D.h:27
Double_t fTsumwz
Definition: TH3.h:37
Histogram is forced to be not weighted even when the histogram is filled with weighted different than...
Definition: TH1.h:155
TArrayD fSumw2
Array of sum of squares of weights.
Definition: TH1.h:94
Double_t * GetB()
Definition: TProfile3D.h:72
virtual TProfile2D * Project3DProfile(Option_t *option="xy") const
Project a 3-D profile into a 2D-profile histogram depending on the option parameter.
Double_t * GetW2()
Definition: TProfile3D.h:75
virtual TProfile2D * Project3DProfile(Option_t *option="xy") const
Project a 3-d histogram into a 2-d profile histograms depending on the option parameter option may co...
Definition: TH3.cxx:2564
virtual Int_t FindGoodLimits(TH1 *h, Double_t xmin, Double_t xmax)
compute the best axis limits for the X axis.
TAxis fZaxis
Z axis descriptor.
Definition: TH1.h:82
Double_t fTsumwyz
Definition: TH3.h:40
Short_t Abs(Short_t d)
Definition: TMathBase.h:108
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 Bool_t Multiply(TF1 *h1, Double_t c1=1)
Performs the operation: this = this*c1*f1 .
Definition: TProfile3D.cxx:855
virtual Double_t GetBinUpEdge(Int_t bin) const
Return up edge of bin.
Definition: TAxis.cxx:514
TH3D()
Constructor.
Definition: TH3.cxx:4174
static void BuildArray(T *p)
virtual TH2D * DoProject2D(const char *name, const char *title, const TAxis *projX, const TAxis *projY, bool computeErrors, bool originalRange, bool useUF, bool useOF) const
internal method performing the projection to a 2D histogram called from TH3::Project3D ...
Definition: TH3.cxx:1899
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
Double_t fTsumwxz
Definition: TH3.h:39
void Class()
Definition: Class.C:29
static Bool_t fgApproximate
Definition: TProfile3D.h:41
virtual void GetStats(Double_t *stats) const
fill the array stats from the contents of this profile.
Definition: TProfile3D.cxx:779
const int ny
Definition: kalman.C:17
static double p2(double t, double a, double b, double c)
virtual Double_t GetBinCenter(Int_t bin) const
Return center of bin.
Definition: TAxis.cxx:464
Double_t fTsumwt
True when TProfile3D::Scale is called.
Definition: TProfile3D.h:38
TString & Append(const char *cs)
Definition: TString.h:497
Double_t * fArray
Definition: TArrayD.h:30
static T * ExtendAxis(T *p, Double_t x, TAxis *axis)
virtual void Copy(TObject &hnew) const
Copy a Profile3D histogram to a new profile2D histogram.
Definition: TProfile3D.cxx:344
virtual Double_t GetBinContent(Int_t bin) const
Return bin content of a Profile3D histogram.
Definition: TProfile3D.cxx:685
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
Double_t fTsumwx
Total Sum of weight*X.
Definition: TH1.h:88
Int_t fN
Definition: TArray.h:38
static void Approximate(Bool_t approx=kTRUE)
Set the fgApproximate flag.
Definition: TProfile3D.cxx:225
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
tomato 3-D histogram with a double per channel (see TH1 documentation)}
Definition: TH3.h:304
Double_t fTsumwy2
Definition: TH3.h:35
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:436
Int_t Fill(const Double_t *v)
Definition: TProfile3D.h:52
Collection abstract base class.
Definition: TCollection.h:42
static Int_t fgBufferSize
!default buffer size for automatic histograms
Definition: TH1.h:104
Double_t fEntries
Number of entries.
Definition: TH1.h:85
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:873
virtual void PutStats(Double_t *stats)
Replace current statistics with the values in array stats.
Definition: TH3.cxx:2647
virtual Bool_t Add(TF1 *h1, Double_t c1=1, Option_t *option="")
Performs the operation: this = this + c1*f1 .
Definition: TProfile3D.cxx:166
TAxis * GetYaxis()
Definition: TH1.h:301
static double p1(double t, double a, double b)
float xmax
Definition: THbookFile.cxx:93
virtual Long64_t Merge(TCollection *list)
Merge all histograms in the collection in this histogram.
Definition: TProfile3D.cxx:847
virtual void SetErrorOption(Option_t *option="")
Set option to compute profile3D errors.
const Bool_t kFALSE
Definition: RtypesCore.h:92
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 Double_t GetBinEffectiveEntries(Int_t bin)
Return bin effective entries for a weighted filled Profile histogram.
Definition: TProfile3D.cxx:715
Double_t fTsumw2
Total Sum of squares of weights.
Definition: TH1.h:87
TProfile3D()
Default constructor for Profile3D histograms.
Definition: TProfile3D.cxx:77
Double_t fTsumwz2
Definition: TH3.h:38
return c2
Definition: legend2.C:14
#define ClassImp(name)
Definition: Rtypes.h:336
virtual void PutStats(Double_t *stats)
Replace current statistics with the values in array stats.
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
Double_t fTsumw
Total Sum of weights.
Definition: TH1.h:86
virtual void SetBinsLength(Int_t n=-1)
Set total number of bins including under/overflow Reallocate bin contents array.
Definition: TH3.cxx:4263
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
TArrayD fBinEntries
Definition: TProfile3D.h:33
void SetBins(const Int_t *nbins, const Double_t *range)
Definition: TProfile3D.h:49
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
Double_t * GetW()
Definition: TProfile3D.h:74
virtual void Copy(TObject &hnew) const
Copy this 3-D histogram structure to newth3.
Definition: TH3.cxx:4242
TAxis * GetZaxis()
Definition: TH1.h:302
virtual void SetBinsLength(Int_t n=-1)
Set total number of bins including under/overflow.
Mother of all ROOT objects.
Definition: TObject.h:37
you should not use this method at all Int_t Int_t z
Definition: TRolke.cxx:630
Bool_t fScaling
Definition: TProfile3D.h:37
Double_t fTsumwxy
Definition: TH3.h:36
EErrorType fErrorMode
Definition: TProfile3D.h:34
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
virtual Int_t BufferEmpty(Int_t action=0)
Fill histogram with all entries in the buffer.
Definition: TProfile3D.cxx:240
virtual void SetEntries(Double_t n)
Definition: TH1.h:363
virtual void ExtendAxis(Double_t x, TAxis *axis)
Profile histogram is resized along axis such that x is in the axis range.
TAxis fXaxis
X axis descriptor.
Definition: TH1.h:80
static void Sumw2(T *p, Bool_t flag)
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 const char * GetName() const
Returns name of object.
Definition: TObject.cxx:364
Double_t * fBuffer
[fBufferSize] entry buffer
Definition: TH1.h:98
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
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)
const TArrayD * GetXbins() const
Definition: TAxis.h:130
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:859
TAxis * GetXaxis()
Definition: TH1.h:300
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
virtual Int_t GetNbinsY() const
Definition: TH1.h:278
virtual Bool_t Divide(TF1 *h1, Double_t c1=1)
Performs the operation: this = this/(c1*f1) .
Definition: TProfile3D.cxx:373
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
virtual TArrayD * GetBinSumw2()
Definition: TProfile2D.h:114
TArrayD fBinSumw2
Definition: TProfile3D.h:40