12 #ifndef ROOT_TProfileHelper
13 #define ROOT_TProfileHelper
30 class TProfileHelper {
34 static Bool_t Add(T* p,
const TH1 *h1,
const TH1 *h2, Double_t c1, Double_t c2=1);
37 static void BuildArray(T* p);
40 static Double_t GetBinEffectiveEntries(T* p, Int_t bin);
43 static Long64_t Merge(T* p, TCollection *list);
46 static T* ExtendAxis(T* p, Double_t x, TAxis *axis);
49 static void Scale(T* p, Double_t c1, Option_t * option);
52 static void Sumw2(T* p, Bool_t flag );
55 static void LabelsDeflate(T* p, Option_t *);
58 static void LabelsInflate(T* p, Option_t *);
61 static Double_t GetBinError(T* p, Int_t bin);
64 static void SetBinEntries(T* p, Int_t bin, Double_t w);
67 static void SetErrorOption(T* p, Option_t * opt);
71 Bool_t TProfileHelper::Add(T* p,
const TH1 *h1,
const TH1 *h2, Double_t c1, Double_t c2)
79 if (p->fBuffer) p->BufferEmpty(1);
82 Int_t nx = p->GetNbinsX();
83 Int_t ny = p->GetNbinsY();
84 Int_t nz = p->GetNbinsZ();
86 if ( nx != p1->GetNbinsX() || nx != p2->GetNbinsX() ||
87 ny != p1->GetNbinsY() || ny != p2->GetNbinsY() ||
88 nz != p1->GetNbinsZ() || nz != p2->GetNbinsZ() ) {
89 Error(
"TProfileHelper::Add",
"Attempt to add profiles with different number of bins");
94 Double_t ac1 = TMath::Abs(c1);
95 Double_t ac2 = TMath::Abs(c2);
96 p->fEntries = ac1*p1->GetEntries() + ac2*p2->GetEntries();
97 Double_t s0[TH1::kNstat], s1[TH1::kNstat], s2[TH1::kNstat];
99 for (i=0;i<TH1::kNstat;i++) {s0[i] = s1[i] = s2[i] = 0;}
103 for (i=0;i<TH1::kNstat;i++) {
104 if (i == 1) s0[i] = c1*c1*s1[i] + c2*c2*s2[i];
105 else s0[i] = ac1*s1[i] + ac2*s2[i];
111 Double_t *cu1 = p1->GetW(); Double_t *cu2 = p2->GetW();
112 Double_t *er1 = p1->GetW2(); Double_t *er2 = p2->GetW2();
113 Double_t *en1 = p1->GetB(); Double_t *en2 = p2->GetB();
114 Double_t *ew1 = p1->GetB2(); Double_t *ew2 = p2->GetB2();
116 if (p->fBinSumw2.fN == 0 && (p1->fBinSumw2.fN != 0 || p2->fBinSumw2.fN != 0) ) p->Sumw2();
118 if (ew1 == 0) ew1 = en1;
119 if (ew2 == 0) ew2 = en2;
120 for (bin =0;bin< p->fN;bin++) {
121 p->fArray[bin] = c1*cu1[bin] + c2*cu2[bin];
122 p->fSumw2.fArray[bin] = ac1*er1[bin] + ac2*er2[bin];
123 p->fBinEntries.fArray[bin] = ac1*en1[bin] + ac2*en2[bin];
124 if (p->fBinSumw2.fN ) p->fBinSumw2.fArray[bin] = ac1*ac1*ew1[bin] + ac2*ac2*ew2[bin];
129 template <
typename T>
130 void TProfileHelper::BuildArray(T* p) {
136 p->fBinEntries.Set(p->fNcells);
137 p->fSumw2.Set(p->fNcells);
138 if (TH1::GetDefaultSumw2() || p->fBinSumw2.fN > 0 ) p->fBinSumw2.Set(p->fNcells);
142 template <
typename T>
143 Double_t TProfileHelper::GetBinEffectiveEntries(T* p, Int_t bin)
153 if (p->fBuffer) p->BufferEmpty();
155 if (bin < 0 || bin >= p->fNcells)
return 0;
156 double sumOfWeights = p->fBinEntries.fArray[bin];
157 if ( p->fBinSumw2.fN == 0 || p->fBinSumw2.fN != p->fNcells) {
162 double sumOfWeightsSquare = p->fBinSumw2.fArray[bin];
163 return ( sumOfWeightsSquare > 0 ? sumOfWeights * sumOfWeights / sumOfWeightsSquare : 0 );
166 template <
typename T>
167 Long64_t TProfileHelper::Merge(T* p, TCollection *li) {
182 if (li->IsEmpty())
return (Int_t) p->GetEntries();
190 Bool_t initialLimitsFound = kFALSE;
191 Bool_t allSameLimits = kTRUE;
192 Bool_t allHaveLimits = kTRUE;
203 Bool_t hasLimits = h->GetXaxis()->GetXmin() < h->GetXaxis()->GetXmax();
204 allHaveLimits = allHaveLimits && hasLimits;
212 if (firstNonEmptyHist ) {
215 if (!p->SameLimitsAndNBins(p->fXaxis, *(h->GetXaxis())) )
216 p->fXaxis.Set(h->GetXaxis()->GetNbins(), h->GetXaxis()->GetXmin(),h->GetXaxis()->GetXmax());
217 if (!p->SameLimitsAndNBins(p->fYaxis, *(h->GetYaxis())) )
218 p->fYaxis.Set(h->GetYaxis()->GetNbins(), h->GetYaxis()->GetXmin(),h->GetYaxis()->GetXmax());
219 if (!p->SameLimitsAndNBins(p->fZaxis, *(h->GetZaxis())) )
220 p->fZaxis.Set(h->GetZaxis()->GetNbins(), h->GetZaxis()->GetXmin(),h->GetZaxis()->GetXmax());
222 firstNonEmptyHist = kFALSE;
228 if (!initialLimitsFound) {
229 initialLimitsFound = kTRUE;
230 newXAxis.Set(h->GetXaxis()->GetNbins(), h->GetXaxis()->GetXmin(),
231 h->GetXaxis()->GetXmax());
232 if ( p->GetDimension() >= 2 )
233 newYAxis.Set(h->GetYaxis()->GetNbins(), h->GetYaxis()->GetXmin(),
234 h->GetYaxis()->GetXmax());
235 if ( p->GetDimension() >= 3 )
236 newZAxis.Set(h->GetZaxis()->GetNbins(), h->GetZaxis()->GetXmin(),
237 h->GetZaxis()->GetXmax());
241 if (!p->SameLimitsAndNBins(newXAxis, *(h->GetXaxis())) ||
242 !p->SameLimitsAndNBins(newYAxis, *(h->GetYaxis())) ||
243 !p->SameLimitsAndNBins(newZAxis, *(h->GetZaxis())) ) {
245 allSameLimits = kFALSE;
250 if (!p->RecomputeAxisLimits(newXAxis, *(h->GetXaxis()))) {
251 Error(
"TProfileHelper::Merge",
"Cannot merge profiles %d dim - limits are inconsistent:\n "
252 "first: (%d, %f, %f), second: (%d, %f, %f)",p->GetDimension(),
253 newXAxis.GetNbins(), newXAxis.GetXmin(), newXAxis.GetXmax(),
254 h->GetXaxis()->GetNbins(), h->GetXaxis()->GetXmin(),
255 h->GetXaxis()->GetXmax());
258 if (p->GetDimension() >= 2 && !p->RecomputeAxisLimits(newYAxis, *(h->GetYaxis()))) {
259 Error(
"TProfileHelper::Merge",
"Cannot merge profiles %d dim - limits are inconsistent:\n "
260 "first: (%d, %f, %f), second: (%d, %f, %f)",p->GetDimension(),
261 newYAxis.GetNbins(), newYAxis.GetXmin(), newYAxis.GetXmax(),
262 h->GetYaxis()->GetNbins(), h->GetYaxis()->GetXmin(),
263 h->GetYaxis()->GetXmax());
266 if (p->GetDimension() >= 3 && !p->RecomputeAxisLimits(newZAxis, *(h->GetZaxis()))) {
267 Error(
"TProfileHelper::Merge",
"Cannot merge profiles %d dim - limits are inconsistent:\n "
268 "first: (%d, %f, %f), second: (%d, %f, %f)",p->GetDimension(),
269 newZAxis.GetNbins(), newZAxis.GetXmin(), newZAxis.GetXmax(),
270 h->GetZaxis()->GetNbins(), h->GetZaxis()->GetXmin(),
271 h->GetZaxis()->GetXmax());
277 }
while ( ( h = dynamic_cast<T*> ( next() ) ) != NULL );
278 if (!h && (*next) ) {
279 Error(
"TProfileHelper::Merge",
"Attempt to merge object of class: %s to a %s",
280 (*next)->ClassName(),p->ClassName());
291 if (!allSameLimits) {
294 Bool_t mustCleanup = p->TestBit(kMustCleanup);
295 if (mustCleanup) p->ResetBit(kMustCleanup);
296 hclone = (T*)p->IsA()->New();
298 hclone->SetDirectory(0);
300 if (mustCleanup) p->SetBit(kMustCleanup);
304 inlist.AddFirst(hclone);
307 if (!allSameLimits && initialLimitsFound) {
308 Int_t b[] = { newXAxis.GetNbins(), newYAxis.GetNbins(), newZAxis.GetNbins() };
309 Double_t v[] = { newXAxis.GetXmin(), newXAxis.GetXmax(),
310 newYAxis.GetXmin(), newYAxis.GetXmax(),
311 newZAxis.GetXmin(), newZAxis.GetXmax() };
315 if (!allHaveLimits) {
317 while ( (h = dynamic_cast<T*> (next()) ) ) {
318 if (h->GetXaxis()->GetXmin() >= h->GetXaxis()->GetXmax() && h->fBuffer) {
320 Int_t nbentries = (Int_t)h->fBuffer[0];
322 for (Int_t i = 0; i < nbentries; i++)
323 if ( p->GetDimension() == 3 ) {
324 v[0] = h->fBuffer[5*i + 2];
325 v[1] = h->fBuffer[5*i + 3];
326 v[2] = h->fBuffer[5*i + 4];
327 v[3] = h->fBuffer[5*i + 5];
328 v[4] = h->fBuffer[5*i + 1];
330 }
else if ( p->GetDimension() == 2 ) {
331 v[0] = h->fBuffer[4*i + 2];
332 v[1] = h->fBuffer[4*i + 3];
333 v[2] = h->fBuffer[4*i + 4];
334 v[3] = h->fBuffer[4*i + 1];
338 else if ( p->GetDimension() == 1 ) {
339 v[0] = h->fBuffer[3*i + 2];
340 v[1] = h->fBuffer[3*i + 3];
341 v[2] = h->fBuffer[3*i + 1];
347 if (!initialLimitsFound) {
349 inlist.Remove(hclone);
352 return (Int_t) p->GetEntries();
358 Double_t stats[TH1::kNstat], totstats[TH1::kNstat];
359 for (Int_t i=0;i<TH1::kNstat;i++) {totstats[i] = stats[i] = 0;}
360 p->GetStats(totstats);
361 Double_t nentries = p->GetEntries();
362 Bool_t canExtend = p->CanExtendAllAxes();
363 p->SetCanExtend(TH1::kNoAxis);
365 while ( (h=static_cast<T*>(next())) ) {
368 if (h->GetXaxis()->GetXmin() < h->GetXaxis()->GetXmax()) {
371 for (Int_t i = 0; i < TH1::kNstat; i++)
372 totstats[i] += stats[i];
373 nentries += h->GetEntries();
375 for ( Int_t hbin = 0; hbin < h->fN; ++hbin ) {
377 if (!allSameLimits) {
382 if ( h->GetW()[hbin] != 0 && (h->IsBinUnderflow(hbin) || h->IsBinOverflow(hbin)) ) {
384 Error(
"TProfileHelper::Merge",
"Cannot merge profiles - they have"
385 " different limits and underflows/overflows are present."
386 " The initial profile is now broken!");
389 Int_t hbinx, hbiny, hbinz;
390 h->GetBinXYZ(hbin, hbinx, hbiny, hbinz);
392 pbin = p->GetBin( p->fXaxis.FindBin( h->GetXaxis()->GetBinCenter(hbinx) ),
393 p->fYaxis.FindBin( h->GetYaxis()->GetBinCenter(hbiny) ),
394 p->fZaxis.FindBin( h->GetZaxis()->GetBinCenter(hbinz) ) );
398 p->fArray[pbin] += h->GetW()[hbin];
399 p->fSumw2.fArray[pbin] += h->GetW2()[hbin];
400 p->fBinEntries.fArray[pbin] += h->GetB()[hbin];
401 if (p->fBinSumw2.fN) {
402 if ( h->GetB2() ) p->fBinSumw2.fArray[pbin] += h->GetB2()[hbin];
403 else p->fBinSumw2.fArray[pbin] += h->GetB()[hbin];
408 if (canExtend) p->SetCanExtend(TH1::kAllAxes);
411 p->PutStats(totstats);
412 p->SetEntries(nentries);
414 inlist.Remove(hclone);
417 return (Long64_t)nentries;
420 template <
typename T>
421 T* TProfileHelper::ExtendAxis(T* p, Double_t x, TAxis *axis)
433 if (!axis->CanExtend())
return 0;
434 if (axis->GetXmin() >= axis->GetXmax())
return 0;
435 if (axis->GetNbins() <= 0)
return 0;
438 if (!p->FindNewAxisLimits(axis, x, xmin, xmax))
442 T* hold = (T*)p->IsA()->New();
444 hold->SetDirectory(0);
447 axis->SetLimits(xmin,xmax);
448 if (p->fBinSumw2.fN) hold->Sumw2();
450 Int_t nbinsx = p->fXaxis.GetNbins();
451 Int_t nbinsy = p->fYaxis.GetNbins();
452 Int_t nbinsz = p->fZaxis.GetNbins();
458 Int_t ix, iy, iz, binx, biny, binz;
459 for (binz=1;binz<=nbinsz;binz++) {
460 bz = hold->GetZaxis()->GetBinCenter(binz);
461 iz = p->fZaxis.FindFixBin(bz);
462 for (biny=1;biny<=nbinsy;biny++) {
463 by = hold->GetYaxis()->GetBinCenter(biny);
464 iy = p->fYaxis.FindFixBin(by);
465 for (binx=1;binx<=nbinsx;binx++) {
466 bx = hold->GetXaxis()->GetBinCenter(binx);
467 ix = p->fXaxis.FindFixBin(bx);
469 Int_t sourceBin = hold->GetBin(binx,biny,binz);
470 Int_t destinationBin = p->GetBin(ix,iy,iz);
471 p->AddBinContent(destinationBin, hold->fArray[sourceBin]);
472 p->fBinEntries.fArray[destinationBin] += hold->fBinEntries.fArray[sourceBin];
473 p->fSumw2.fArray[destinationBin] += hold->fSumw2.fArray[sourceBin];
474 if (p->fBinSumw2.fN) p->fBinSumw2.fArray[destinationBin] += hold->fBinSumw2.fArray[sourceBin];
481 template <
typename T>
482 void TProfileHelper::Scale(T* p, Double_t c1, Option_t *)
484 Double_t ac1 = TMath::Abs(c1);
488 Double_t *cu1 = p->GetW();
489 Double_t *er1 = p->GetW2();
490 Double_t *en1 = p->GetB();
491 for (bin=0;bin<p->fN;bin++) {
492 p->fArray[bin] = c1*cu1[bin];
493 p->fSumw2.fArray[bin] = ac1*ac1*er1[bin];
494 p->fBinEntries.fArray[bin] = en1[bin];
498 template <
typename T>
499 void TProfileHelper::Sumw2(T* p, Bool_t flag)
511 if (p->fBinSumw2.fN > 0 ) p->fBinSumw2.Set(0);
515 if ( p->fBinSumw2.fN == p->fNcells) {
516 if (!p->fgDefaultSumw2)
517 Warning(
"Sumw2",
"Sum of squares of profile bin weights structure already created");
521 p->fBinSumw2.Set(p->fNcells);
524 for (Int_t bin=0; bin<p->fNcells; bin++) {
525 p->fBinSumw2.fArray[bin] = p->fBinEntries.fArray[bin];
529 template <
typename T>
530 void TProfileHelper::LabelsDeflate(T* p, Option_t *ax)
538 TAxis *axis = p->GetXaxis();
539 if (ax[0] ==
'y' || ax[0] ==
'Y') axis = p->GetYaxis();
540 if (ax[0] ==
'z' || ax[0] ==
'Z') axis = p->GetZaxis();
542 Error(
"TProfileHelper::LabelsDeflate",
"Invalid axis option %s",ax);
545 if (!axis->GetLabels())
return;
550 TIter next(axis->GetLabels());
553 while ((obj = next())) {
554 Int_t ibin = obj->GetUniqueID();
555 if (ibin > nbins) nbins = ibin;
557 if (nbins < 1) nbins = 1;
560 if (nbins==axis->GetNbins())
return;
562 T *hold = (T*)p->IsA()->New();;
563 hold->SetDirectory(0);
567 Double_t xmin = axis->GetXmin();
568 Double_t xmax = axis->GetBinUpEdge(nbins);
571 axis->Set(nbins,xmin,xmax);
572 p->SetBinsLength(-1);
573 p->fBinEntries.Set(p->fN);
574 p->fSumw2.Set(p->fN);
575 if (p->fBinSumw2.fN) p->fBinSumw2.Set(p->fN);
581 Int_t bin,binx,biny,binz;
582 for (bin =0; bin < hold->fN; ++bin)
584 hold->GetBinXYZ(bin, binx, biny, binz);
585 Int_t ibin = p->GetBin(binx, biny, binz);
586 p->fArray[ibin] += hold->fArray[bin];
587 p->fBinEntries.fArray[ibin] += hold->fBinEntries.fArray[bin];
588 p->fSumw2.fArray[ibin] += hold->fSumw2.fArray[bin];
589 if (p->fBinSumw2.fN) p->fBinSumw2.fArray[ibin] += hold->fBinSumw2.fArray[bin];
595 template <
typename T>
596 void TProfileHelper::LabelsInflate(T* p, Option_t *ax)
604 TAxis *axis = p->GetXaxis();
605 if (ax[0] ==
'y' || ax[0] ==
'Y') axis = p->GetYaxis();
606 T *hold = (T*)p->IsA()->New();;
607 hold->SetDirectory(0);
613 Int_t nbins = axis->GetNbins();
614 Double_t xmin = axis->GetXmin();
615 Double_t xmax = axis->GetXmax();
616 xmax = xmin + 2*(xmax-xmin);
620 axis->Set(nbins,xmin,xmax);
622 p->SetBinsLength(-1);
623 Int_t ncells = p->fN;
624 p->fBinEntries.Set(ncells);
625 p->fSumw2.Set(ncells);
626 if (p->fBinSumw2.fN) p->fBinSumw2.Set(ncells);
629 for (Int_t ibin =0; ibin < p->fN; ibin++)
631 Int_t binx, biny, binz;
632 p->GetBinXYZ(ibin, binx, biny, binz);
633 Int_t bin = hold->GetBin(binx, biny, binz);
635 if (p->IsBinUnderflow(ibin) || p->IsBinOverflow(ibin)) {
636 p->UpdateBinContent(ibin, 0.0);
637 p->fBinEntries.fArray[ibin] = 0.0;
638 p->fSumw2.fArray[ibin] = 0.0;
639 if (p->fBinSumw2.fN) p->fBinSumw2.fArray[ibin] = 0.0;
641 p->fArray[ibin] = hold->fArray[bin];
642 p->fBinEntries.fArray[ibin] = hold->fBinEntries.fArray[bin];
643 p->fSumw2.fArray[ibin] = hold->fSumw2.fArray[bin];
644 if (p->fBinSumw2.fN) p->fBinSumw2.fArray[ibin] = hold->fBinSumw2.fArray[bin];
650 template <
typename T>
651 void TProfileHelper::SetErrorOption(T* p, Option_t * option) {
653 TString opt = option;
655 p->fErrorMode = kERRORMEAN;
656 if (opt.Contains(
"s")) p->fErrorMode = kERRORSPREAD;
657 if (opt.Contains(
"i")) p->fErrorMode = kERRORSPREADI;
658 if (opt.Contains(
"g")) p->fErrorMode = kERRORSPREADG;
661 template <
typename T>
662 Double_t TProfileHelper::GetBinError(T* p, Int_t bin)
666 if (p->fBuffer) p->BufferEmpty();
668 if (bin < 0 || bin >= p->fNcells)
return 0;
669 Double_t cont = p->fArray[bin];
670 Double_t sum = p->fBinEntries.fArray[bin];
671 Double_t err2 = p->fSumw2.fArray[bin];
672 Double_t neff = p->GetBinEffectiveEntries(bin);
673 if (sum == 0)
return 0;
675 if (p->fErrorMode == kERRORSPREADG) {
676 return 1./TMath::Sqrt(sum);
679 Double_t contsum = cont/sum;
680 Double_t eprim2 = TMath::Abs(err2/sum - contsum*contsum);
681 Double_t eprim = TMath::Sqrt(eprim2);
683 if (p->fErrorMode == kERRORSPREADI) {
684 if (eprim != 0)
return eprim/TMath::Sqrt(neff);
687 return 1/TMath::Sqrt(12*neff);
693 if (err2 != 0 && neff < 5) test = eprim2*sum/err2;
695 if (p->fgApproximate && (test < 1.e-4 || eprim2 <= 0)) {
696 Double_t stats[TH1::kNstat];
698 Double_t ssum = stats[0];
701 if (p->GetDimension() == 2) index = 7;
702 if (p->GetDimension() == 3) index = 11;
703 Double_t scont = stats[index];
704 Double_t serr2 = stats[index+1];
707 Double_t scontsum = scont/ssum;
708 Double_t seprim2 = TMath::Abs(serr2/ssum - scontsum*scontsum);
709 eprim = 2*TMath::Sqrt(seprim2);
712 sum = TMath::Abs(sum);
715 if (p->fErrorMode == kERRORSPREAD)
return eprim;
720 return eprim/TMath::Sqrt(neff);
725 template <
typename T>
726 void TProfileHelper::SetBinEntries(T* p, Int_t bin, Double_t w) {
731 if (bin < 0 || bin >= p->fNcells)
return;
732 p->fBinEntries.fArray[bin] = w;
733 if (p->fBinSumw2.fN) p->fBinSumw2.fArray[bin] = w;