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