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