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