64 fTsumwy = fTsumwy2 = fTsumwxy = 0;
71 TH2::TH2(
const char *name,
const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
72 ,Int_t nbinsy,Double_t ylow,Double_t yup)
73 :TH1(name,title,nbinsx,xlow,xup)
77 fTsumwy = fTsumwy2 = fTsumwxy = 0;
78 if (nbinsy <= 0) {Warning(
"TH2",
"nbinsy is <=0 - set to nbinsy = 1"); nbinsy = 1; }
79 fYaxis.Set(nbinsy,ylow,yup);
80 fNcells = fNcells*(nbinsy+2);
87 TH2::TH2(
const char *name,
const char *title,Int_t nbinsx,
const Double_t *xbins
88 ,Int_t nbinsy,Double_t ylow,Double_t yup)
89 :TH1(name,title,nbinsx,xbins)
93 fTsumwy = fTsumwy2 = fTsumwxy = 0;
94 if (nbinsy <= 0) {Warning(
"TH2",
"nbinsy is <=0 - set to nbinsy = 1"); nbinsy = 1; }
95 fYaxis.Set(nbinsy,ylow,yup);
96 fNcells = fNcells*(nbinsy+2);
103 TH2::TH2(
const char *name,
const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
104 ,Int_t nbinsy,
const Double_t *ybins)
105 :TH1(name,title,nbinsx,xlow,xup)
109 fTsumwy = fTsumwy2 = fTsumwxy = 0;
110 if (nbinsy <= 0) {Warning(
"TH2",
"nbinsy is <=0 - set to nbinsy = 1"); nbinsy = 1; }
111 if (ybins) fYaxis.Set(nbinsy,ybins);
112 else fYaxis.Set(nbinsy,0,1);
113 fNcells = fNcells*(nbinsy+2);
120 TH2::TH2(
const char *name,
const char *title,Int_t nbinsx,
const Double_t *xbins
121 ,Int_t nbinsy,
const Double_t *ybins)
122 :TH1(name,title,nbinsx,xbins)
126 fTsumwy = fTsumwy2 = fTsumwxy = 0;
127 if (nbinsy <= 0) {Warning(
"TH2",
"nbinsy is <=0 - set to nbinsy = 1"); nbinsy = 1; }
128 if (ybins) fYaxis.Set(nbinsy,ybins);
129 else fYaxis.Set(nbinsy,0,1);
130 fNcells = fNcells*(nbinsy+2);
137 TH2::TH2(
const char *name,
const char *title,Int_t nbinsx,
const Float_t *xbins
138 ,Int_t nbinsy,
const Float_t *ybins)
139 :TH1(name,title,nbinsx,xbins)
143 fTsumwy = fTsumwy2 = fTsumwxy = 0;
144 if (nbinsy <= 0) {Warning(
"TH2",
"nbinsy is <=0 - set to nbinsy = 1"); nbinsy = 1; }
145 if (ybins) fYaxis.Set(nbinsy,ybins);
146 else fYaxis.Set(nbinsy,0,1);
147 fNcells = fNcells*(nbinsy+2);
155 TH2::TH2(
const TH2 &h) : TH1()
157 ((TH2&)h).Copy(*
this);
177 Int_t TH2::BufferEmpty(Int_t action)
180 if (!fBuffer)
return 0;
181 Int_t nbentries = (Int_t)fBuffer[0];
185 if (nbentries == 0)
return 0;
186 if (nbentries < 0 && action == 0)
return 0;
188 Double_t *buffer = fBuffer;
190 nbentries = -nbentries;
199 if (CanExtendAllAxes() || fXaxis.GetXmax() <= fXaxis.GetXmin() || fYaxis.GetXmax() <= fYaxis.GetXmin()) {
201 Double_t xmin = fBuffer[2];
202 Double_t xmax = xmin;
203 Double_t ymin = fBuffer[3];
204 Double_t ymax = ymin;
205 for (Int_t i=1;i<nbentries;i++) {
206 Double_t x = fBuffer[3*i+2];
207 if (x < xmin) xmin = x;
208 if (x > xmax) xmax = x;
209 Double_t y = fBuffer[3*i+3];
210 if (y < ymin) ymin = y;
211 if (y > ymax) ymax = y;
213 if (fXaxis.GetXmax() <= fXaxis.GetXmin() || fYaxis.GetXmax() <= fYaxis.GetXmin()) {
214 THLimitsFinder::GetLimitsFinder()->FindGoodLimits(
this,xmin,xmax,ymin,ymax);
217 Int_t keep = fBufferSize; fBufferSize = 0;
218 if (xmin < fXaxis.GetXmin()) ExtendAxis(xmin,&fXaxis);
219 if (xmax >= fXaxis.GetXmax()) ExtendAxis(xmax,&fXaxis);
220 if (ymin < fYaxis.GetXmin()) ExtendAxis(ymin,&fYaxis);
221 if (ymax >= fYaxis.GetXmax()) ExtendAxis(ymax,&fYaxis);
228 for (Int_t i=0;i<nbentries;i++) {
229 Fill(buffer[3*i+2],buffer[3*i+3],buffer[3*i+1]);
233 if (action > 0) {
delete [] fBuffer; fBuffer = 0; fBufferSize = 0;}
235 if (nbentries == (Int_t)fEntries) fBuffer[0] = -nbentries;
251 Int_t TH2::BufferFill(Double_t x, Double_t y, Double_t w)
253 if (!fBuffer)
return -3;
254 Int_t nbentries = (Int_t)fBuffer[0];
256 nbentries = -nbentries;
257 fBuffer[0] = nbentries;
259 Double_t *buffer = fBuffer; fBuffer=0;
264 if (3*nbentries+3 >= fBufferSize) {
268 fBuffer[3*nbentries+1] = w;
269 fBuffer[3*nbentries+2] = x;
270 fBuffer[3*nbentries+3] = y;
279 void TH2::Copy(TObject &obj)
const
282 ((TH2&)obj).fScalefactor = fScalefactor;
283 ((TH2&)obj).fTsumwy = fTsumwy;
284 ((TH2&)obj).fTsumwy2 = fTsumwy2;
285 ((TH2&)obj).fTsumwxy = fTsumwxy;
292 Int_t TH2::Fill(Double_t )
294 Error(
"Fill",
"Invalid signature - do nothing");
314 Int_t TH2::Fill(Double_t x,Double_t y)
316 if (fBuffer)
return BufferFill(x,y,1);
318 Int_t binx, biny, bin;
320 binx = fXaxis.FindBin(x);
321 biny = fYaxis.FindBin(y);
322 if (binx <0 || biny <0)
return -1;
323 bin = biny*(fXaxis.GetNbins()+2) + binx;
325 if (fSumw2.fN) ++fSumw2.fArray[bin];
326 if (binx == 0 || binx > fXaxis.GetNbins()) {
327 if (!GetStatOverflowsBehaviour())
return -1;
329 if (biny == 0 || biny > fYaxis.GetNbins()) {
330 if (!GetStatOverflowsBehaviour())
return -1;
358 Int_t TH2::Fill(Double_t x, Double_t y, Double_t w)
360 if (fBuffer)
return BufferFill(x,y,w);
362 Int_t binx, biny, bin;
364 binx = fXaxis.FindBin(x);
365 biny = fYaxis.FindBin(y);
366 if (binx <0 || biny <0)
return -1;
367 bin = biny*(fXaxis.GetNbins()+2) + binx;
368 if (!fSumw2.fN && w != 1.0 && !TestBit(TH1::kIsNotW)) Sumw2();
369 if (fSumw2.fN) fSumw2.fArray[bin] += w*w;
370 AddBinContent(bin,w);
371 if (binx == 0 || binx > fXaxis.GetNbins()) {
372 if (!GetStatOverflowsBehaviour())
return -1;
374 if (biny == 0 || biny > fYaxis.GetNbins()) {
375 if (!GetStatOverflowsBehaviour())
return -1;
404 Int_t TH2::Fill(
const char *namex,
const char *namey, Double_t w)
406 Int_t binx, biny, bin;
408 binx = fXaxis.FindBin(namex);
409 biny = fYaxis.FindBin(namey);
410 if (binx <0 || biny <0)
return -1;
411 bin = biny*(fXaxis.GetNbins()+2) + binx;
412 if (!fSumw2.fN && w != 1.0 && !TestBit(TH1::kIsNotW)) Sumw2();
413 if (fSumw2.fN) fSumw2.fArray[bin] += w*w;
414 AddBinContent(bin,w);
415 if (binx == 0 || binx > fXaxis.GetNbins())
return -1;
416 if (biny == 0 || biny > fYaxis.GetNbins())
return -1;
417 Double_t x = fXaxis.GetBinCenter(binx);
418 Double_t y = fYaxis.GetBinCenter(biny);
446 Int_t TH2::Fill(
const char *namex, Double_t y, Double_t w)
448 Int_t binx, biny, bin;
450 binx = fXaxis.FindBin(namex);
451 biny = fYaxis.FindBin(y);
452 if (binx <0 || biny <0)
return -1;
453 bin = biny*(fXaxis.GetNbins()+2) + binx;
454 if (!fSumw2.fN && w != 1.0 && !TestBit(TH1::kIsNotW)) Sumw2();
455 if (fSumw2.fN) fSumw2.fArray[bin] += w*w;
456 AddBinContent(bin,w);
457 if (binx == 0 || binx > fXaxis.GetNbins())
return -1;
458 if (biny == 0 || biny > fYaxis.GetNbins()) {
459 if (!GetStatOverflowsBehaviour())
return -1;
461 Double_t x = fXaxis.GetBinCenter(binx);
489 Int_t TH2::Fill(Double_t x,
const char *namey, Double_t w)
491 Int_t binx, biny, bin;
493 binx = fXaxis.FindBin(x);
494 biny = fYaxis.FindBin(namey);
495 if (binx <0 || biny <0)
return -1;
496 bin = biny*(fXaxis.GetNbins()+2) + binx;
497 if (!fSumw2.fN && w != 1.0 && !TestBit(TH1::kIsNotW)) Sumw2();
498 if (fSumw2.fN) fSumw2.fArray[bin] += w*w;
499 AddBinContent(bin,w);
500 if (binx == 0 || binx > fXaxis.GetNbins()) {
501 if (!GetStatOverflowsBehaviour())
return -1;
503 if (biny == 0 || biny > fYaxis.GetNbins())
return -1;
504 Double_t y = fYaxis.GetBinCenter(biny);
533 void TH2::FillN(Int_t ntimes,
const Double_t *x,
const Double_t *y,
const Double_t *w, Int_t stride)
535 Int_t binx, biny, bin, i;
542 for (i=0;i<ntimes;i+=stride) {
544 if (w) BufferFill(x[i],y[i],w[i]);
545 else BufferFill(x[i], y[i], 1.);
548 if (i < ntimes && fBuffer==0)
555 for (i=ifirst;i<ntimes;i+=stride) {
557 binx = fXaxis.FindBin(x[i]);
558 biny = fYaxis.FindBin(y[i]);
559 if (binx <0 || biny <0)
continue;
560 bin = biny*(fXaxis.GetNbins()+2) + binx;
562 if (!fSumw2.fN && ww != 1.0 && !TestBit(TH1::kIsNotW)) Sumw2();
563 if (fSumw2.fN) fSumw2.fArray[bin] += ww*ww;
564 AddBinContent(bin,ww);
565 if (binx == 0 || binx > fXaxis.GetNbins()) {
566 if (!GetStatOverflowsBehaviour())
continue;
568 if (biny == 0 || biny > fYaxis.GetNbins()) {
569 if (!GetStatOverflowsBehaviour())
continue;
575 fTsumwx2 += z*x[i]*x[i];
577 fTsumwy2 += z*y[i]*y[i];
578 fTsumwxy += z*x[i]*y[i];
597 void TH2::FillRandom(
const char *fname, Int_t ntimes)
599 Int_t bin, binx, biny, ibin, loop;
602 TObject *fobj = gROOT->GetFunction(fname);
603 if (!fobj) { Error(
"FillRandom",
"Unknown function: %s",fname);
return; }
604 TF2 * f1 =
dynamic_cast<TF2*
>(fobj);
605 if (!f1) { Error(
"FillRandom",
"Function: %s is not a TF2, is a %s",fname,fobj->IsA()->GetName());
return; }
608 TAxis & xAxis = fXaxis;
609 TAxis & yAxis = fYaxis;
612 if (fXaxis.GetXmax() <= fXaxis.GetXmin() || fYaxis.GetXmax() <= fYaxis.GetXmin()) {
613 Double_t xmin,xmax,ymin,ymax;
614 f1->GetRange(xmin,ymin,xmax,ymax);
615 Info(
"FillRandom",
"Using function axis and range ([%g,%g],[%g,%g])",xmin, xmax,ymin,ymax);
616 xAxis = *(f1->GetHistogram()->GetXaxis());
617 yAxis = *(f1->GetHistogram()->GetYaxis());
622 Int_t nbinsx = xAxis.GetNbins();
623 Int_t nbinsy = yAxis.GetNbins();
624 Int_t nbins = nbinsx*nbinsy;
627 Double_t *integral =
new Double_t[nbins+1];
630 for (biny=1;biny<=nbinsy;biny++) {
631 for (binx=1;binx<=nbinsx;binx++) {
633 Double_t fint = f1->Integral(xAxis.GetBinLowEdge(binx), xAxis.GetBinUpEdge(binx), yAxis.GetBinLowEdge(biny), yAxis.GetBinUpEdge(biny));
634 integral[ibin] = integral[ibin-1] + fint;
639 if (integral[nbins] == 0 ) {
641 Error(
"FillRandom",
"Integral = zero");
return;
643 for (bin=1;bin<=nbins;bin++) integral[bin] /= integral[nbins];
646 for (loop=0;loop<ntimes;loop++) {
647 r1 = gRandom->Rndm();
648 ibin = TMath::BinarySearch(nbins,&integral[0],r1);
650 binx = 1 + ibin - nbinsx*biny;
652 x = xAxis.GetBinCenter(binx);
653 y = yAxis.GetBinCenter(biny);
672 void TH2::FillRandom(TH1 *h, Int_t ntimes)
674 if (!h) { Error(
"FillRandom",
"Null histogram");
return; }
675 if (fDimension != h->GetDimension()) {
676 Error(
"FillRandom",
"Histograms with different dimensions");
return;
679 if (h->ComputeIntegral() == 0)
return;
684 for (loop=0;loop<ntimes;loop++) {
693 void TH2::DoFitSlices(
bool onX,
694 TF1 *f1, Int_t firstbin, Int_t lastbin, Int_t cut, Option_t *option, TObjArray* arr)
696 TAxis& outerAxis = (onX ? fYaxis : fXaxis);
697 TAxis& innerAxis = (onX ? fXaxis : fYaxis);
699 Int_t nbins = outerAxis.GetNbins();
700 if (firstbin < 0) firstbin = 0;
701 if (lastbin < 0 || lastbin > nbins + 1) lastbin = nbins + 1;
702 if (lastbin < firstbin) {firstbin = 0; lastbin = nbins + 1;}
703 TString opt = option;
704 TString proj_opt =
"e";
705 Int_t i1 = opt.Index(
"[");
706 Int_t i2 = opt.Index(
"]");
707 if (i1>=0 && i2>i1) {
708 proj_opt += opt(i1,i2-i1+1);
709 opt.Remove(i1, i2-i1+1);
713 if (opt.Contains(
"g2")) {ngroup = 2; opt.ReplaceAll(
"g2",
"");}
714 if (opt.Contains(
"g3")) {ngroup = 3; opt.ReplaceAll(
"g3",
"");}
715 if (opt.Contains(
"g4")) {ngroup = 4; opt.ReplaceAll(
"g4",
"");}
716 if (opt.Contains(
"g5")) {ngroup = 5; opt.ReplaceAll(
"g5",
"");}
719 Int_t nstep = ngroup;
720 if (opt.Contains(
"s")) nstep = 1;
724 f1 = (TF1*)gROOT->GetFunction(
"gaus");
725 if (f1 == 0) f1 =
new TF1(
"gaus",
"gaus",innerAxis.GetXmin(),innerAxis.GetXmax());
726 else f1->SetRange(innerAxis.GetXmin(),innerAxis.GetXmax());
728 Int_t npar = f1->GetNpar();
729 if (npar <= 0)
return;
730 Double_t *parsave =
new Double_t[npar];
731 f1->GetParameters(parsave);
735 arr->Expand(npar + 1);
740 TH1D **hlist =
new TH1D*[npar];
741 char *name =
new char[2000];
742 char *title =
new char[2000];
743 const TArrayD *bins = outerAxis.GetXbins();
744 for (ipar=0;ipar<npar;ipar++) {
745 snprintf(name,2000,
"%s_%d",GetName(),ipar);
746 snprintf(title,2000,
"Fitted value of par[%d]=%s",ipar,f1->GetParName(ipar));
747 delete gDirectory->FindObject(name);
749 hlist[ipar] =
new TH1D(name,title, nbins, outerAxis.GetXmin(), outerAxis.GetXmax());
751 hlist[ipar] =
new TH1D(name,title, nbins,bins->fArray);
753 hlist[ipar]->GetXaxis()->SetTitle(outerAxis.GetTitle());
755 (*arr)[ipar] = hlist[ipar];
757 snprintf(name,2000,
"%s_chi2",GetName());
758 delete gDirectory->FindObject(name);
761 hchi2 =
new TH1D(name,
"chisquare", nbins, outerAxis.GetXmin(), outerAxis.GetXmax());
763 hchi2 =
new TH1D(name,
"chisquare", nbins, bins->fArray);
765 hchi2->GetXaxis()->SetTitle(outerAxis.GetTitle());
767 (*arr)[npar] = hchi2;
774 for (bin=firstbin;bin+ngroup-1<=lastbin;bin += nstep) {
777 hp= ProjectionX(
"_temp",bin,bin+ngroup-1,proj_opt);
779 hp= ProjectionY(
"_temp",bin,bin+ngroup-1,proj_opt);
780 if (hp == 0)
continue;
781 nentries = Long64_t(hp->GetEntries());
782 if (nentries == 0 || nentries < cut) {
delete hp;
continue;}
783 f1->SetParameters(parsave);
784 hp->Fit(f1,opt.Data());
785 Int_t npfits = f1->GetNumberFitPoints();
786 if (npfits > npar && npfits >= cut) {
787 Int_t binOn = bin + ngroup/2;
788 for (ipar=0;ipar<npar;ipar++) {
789 hlist[ipar]->Fill(outerAxis.GetBinCenter(binOn),f1->GetParameter(ipar));
790 hlist[ipar]->SetBinError(binOn,f1->GetParError(ipar));
792 hchi2->SetBinContent(binOn,f1->GetChisquare()/(npfits-npar));
858 void TH2::FitSlicesX(TF1 *f1, Int_t firstybin, Int_t lastybin, Int_t cut, Option_t *option, TObjArray* arr)
860 DoFitSlices(
true, f1, firstybin, lastybin, cut, option, arr);
923 void TH2::FitSlicesY(TF1 *f1, Int_t firstxbin, Int_t lastxbin, Int_t cut, Option_t *option, TObjArray* arr)
925 DoFitSlices(
false, f1, firstxbin, lastxbin, cut, option, arr);
928 Int_t TH2::GetBin(Int_t binx, Int_t biny, Int_t)
const
931 Int_t ofy = fYaxis.GetNbins() + 1;
932 if (biny < 0) biny = 0;
933 if (biny > ofy) biny = ofy;
935 return TH1::GetBin(binx) + (fXaxis.GetNbins() + 2) * biny;
959 Double_t TH2::GetBinWithContent2(Double_t c, Int_t &binx, Int_t &biny, Int_t firstxbin, Int_t lastxbin,
960 Int_t firstybin, Int_t lastybin, Double_t maxdiff)
const
962 if (fDimension != 2) {
965 Error(
"GetBinWithContent2",
"function is only valid for 2-D histograms");
968 if (firstxbin < 0) firstxbin = 1;
969 if (lastxbin < firstxbin) lastxbin = fXaxis.GetNbins();
970 if (firstybin < 0) firstybin = 1;
971 if (lastybin < firstybin) lastybin = fYaxis.GetNbins();
972 Double_t diff, curmax = 1.e240;
973 for (Int_t j = firstybin; j <= lastybin; j++) {
974 for (Int_t i = firstxbin; i <= lastxbin; i++) {
975 diff = TMath::Abs(GetBinContent(i,j)-c);
976 if (diff <= 0) {binx = i; biny=j;
return diff;}
977 if (diff < curmax && diff <= maxdiff) {curmax = diff, binx=i; biny=j;}
987 Double_t TH2::GetCorrelationFactor(Int_t axis1, Int_t axis2)
const
989 if (axis1 < 1 || axis2 < 1 || axis1 > 2 || axis2 > 2) {
990 Error(
"GetCorrelationFactor",
"Wrong parameters");
993 if (axis1 == axis2)
return 1;
994 Double_t stddev1 = GetStdDev(axis1);
995 if (stddev1 == 0)
return 0;
996 Double_t stddev2 = GetStdDev(axis2);
997 if (stddev2 == 0)
return 0;
998 return GetCovariance(axis1,axis2)/stddev1/stddev2;
1005 Double_t TH2::GetCovariance(Int_t axis1, Int_t axis2)
const
1007 if (axis1 < 1 || axis2 < 1 || axis1 > 2 || axis2 > 2) {
1008 Error(
"GetCovariance",
"Wrong parameters");
1011 Double_t stats[kNstat];
1013 Double_t sumw = stats[0];
1015 Double_t sumwx = stats[2];
1016 Double_t sumwx2 = stats[3];
1017 Double_t sumwy = stats[4];
1018 Double_t sumwy2 = stats[5];
1019 Double_t sumwxy = stats[6];
1021 if (sumw == 0)
return 0;
1022 if (axis1 == 1 && axis2 == 1) {
1023 return TMath::Abs(sumwx2/sumw - sumwx/sumw*sumwx/sumw);
1025 if (axis1 == 2 && axis2 == 2) {
1026 return TMath::Abs(sumwy2/sumw - sumwy/sumw*sumwy/sumw);
1028 return sumwxy/sumw - sumwx/sumw*sumwy/sumw;
1037 void TH2::GetRandom2(Double_t &x, Double_t &y)
1039 Int_t nbinsx = GetNbinsX();
1040 Int_t nbinsy = GetNbinsY();
1041 Int_t nbins = nbinsx*nbinsy;
1045 if (fIntegral[nbins+1] != fEntries) integral = ComputeIntegral(
true);
1046 else integral = fIntegral[nbins];
1048 integral = ComputeIntegral(
true);
1050 if (integral == 0 ) { x = 0; y = 0;
return;}
1052 if (integral == TMath::QuietNaN() ) { x = TMath::QuietNaN(); y = TMath::QuietNaN();
return;}
1054 Double_t r1 = gRandom->Rndm();
1055 Int_t ibin = TMath::BinarySearch(nbins,fIntegral,(Double_t) r1);
1056 Int_t biny = ibin/nbinsx;
1057 Int_t binx = ibin - nbinsx*biny;
1058 x = fXaxis.GetBinLowEdge(binx+1);
1059 if (r1 > fIntegral[ibin]) x +=
1060 fXaxis.GetBinWidth(binx+1)*(r1-fIntegral[ibin])/(fIntegral[ibin+1] - fIntegral[ibin]);
1061 y = fYaxis.GetBinLowEdge(biny+1) + fYaxis.GetBinWidth(biny+1)*gRandom->Rndm();
1090 void TH2::GetStats(Double_t *stats)
const
1092 if (fBuffer) ((TH2*)
this)->BufferEmpty();
1094 if ((fTsumw == 0 && fEntries > 0) || fXaxis.TestBit(TAxis::kAxisRange) || fYaxis.TestBit(TAxis::kAxisRange)) {
1095 std::fill(stats, stats + 7, 0);
1097 Int_t firstBinX = fXaxis.GetFirst();
1098 Int_t lastBinX = fXaxis.GetLast();
1099 Int_t firstBinY = fYaxis.GetFirst();
1100 Int_t lastBinY = fYaxis.GetLast();
1102 if (GetStatOverflowsBehaviour()) {
1103 if ( !fXaxis.TestBit(TAxis::kAxisRange) ) {
1104 if (firstBinX == 1) firstBinX = 0;
1105 if (lastBinX == fXaxis.GetNbins() ) lastBinX += 1;
1107 if ( !fYaxis.TestBit(TAxis::kAxisRange) ) {
1108 if (firstBinY == 1) firstBinY = 0;
1109 if (lastBinY == fYaxis.GetNbins() ) lastBinY += 1;
1112 for (Int_t biny = firstBinY; biny <= lastBinY; ++biny) {
1113 Double_t y = fYaxis.GetBinCenter(biny);
1114 for (Int_t binx = firstBinX; binx <= lastBinX; ++binx) {
1115 Double_t x = fXaxis.GetBinCenter(binx);
1117 Int_t bin = GetBin(binx,biny);
1118 Double_t w = RetrieveBinContent(bin);
1119 Double_t wx = w * x;
1120 Double_t wy = w * y;
1123 stats[1] += GetBinErrorSqUnchecked(bin);
1135 stats[3] = fTsumwx2;
1137 stats[5] = fTsumwy2;
1138 stats[6] = fTsumwxy;
1149 Double_t TH2::Integral(Option_t *option)
const
1151 return Integral(fXaxis.GetFirst(),fXaxis.GetLast(),
1152 fYaxis.GetFirst(),fYaxis.GetLast(),option);
1163 Double_t TH2::Integral(Int_t firstxbin, Int_t lastxbin, Int_t firstybin, Int_t lastybin, Option_t *option)
const
1166 return DoIntegral(firstxbin,lastxbin,firstybin,lastybin,-1,0,err,option);
1177 Double_t TH2::IntegralAndError(Int_t firstxbin, Int_t lastxbin, Int_t firstybin, Int_t lastybin, Double_t & error, Option_t *option)
const
1179 return DoIntegral(firstxbin,lastxbin,firstybin,lastybin,-1,0,error,option,kTRUE);
1185 Double_t TH2::Interpolate(Double_t)
const
1187 Error(
"Interpolate",
"This function must be called with 2 arguments for a TH2");
1198 Double_t TH2::Interpolate(Double_t x, Double_t y)
const
1201 Double_t x1=0,x2=0,y1=0,y2=0;
1203 Int_t bin_x = fXaxis.FindFixBin(x);
1204 Int_t bin_y = fYaxis.FindFixBin(y);
1205 if(bin_x<1 || bin_x>GetNbinsX() || bin_y<1 || bin_y>GetNbinsY()) {
1206 Error(
"Interpolate",
"Cannot interpolate outside histogram domain.");
1211 dx = fXaxis.GetBinUpEdge(bin_x)-x;
1212 dy = fYaxis.GetBinUpEdge(bin_y)-y;
1213 if (dx<=fXaxis.GetBinWidth(bin_x)/2 && dy<=fYaxis.GetBinWidth(bin_y)/2)
1215 if (dx>fXaxis.GetBinWidth(bin_x)/2 && dy<=fYaxis.GetBinWidth(bin_y)/2)
1217 if (dx>fXaxis.GetBinWidth(bin_x)/2 && dy>fYaxis.GetBinWidth(bin_y)/2)
1219 if (dx<=fXaxis.GetBinWidth(bin_x)/2 && dy>fYaxis.GetBinWidth(bin_y)/2)
1223 x1 = fXaxis.GetBinCenter(bin_x);
1224 y1 = fYaxis.GetBinCenter(bin_y);
1225 x2 = fXaxis.GetBinCenter(bin_x+1);
1226 y2 = fYaxis.GetBinCenter(bin_y+1);
1229 x1 = fXaxis.GetBinCenter(bin_x-1);
1230 y1 = fYaxis.GetBinCenter(bin_y);
1231 x2 = fXaxis.GetBinCenter(bin_x);
1232 y2 = fYaxis.GetBinCenter(bin_y+1);
1235 x1 = fXaxis.GetBinCenter(bin_x-1);
1236 y1 = fYaxis.GetBinCenter(bin_y-1);
1237 x2 = fXaxis.GetBinCenter(bin_x);
1238 y2 = fYaxis.GetBinCenter(bin_y);
1241 x1 = fXaxis.GetBinCenter(bin_x);
1242 y1 = fYaxis.GetBinCenter(bin_y-1);
1243 x2 = fXaxis.GetBinCenter(bin_x+1);
1244 y2 = fYaxis.GetBinCenter(bin_y);
1247 Int_t bin_x1 = fXaxis.FindFixBin(x1);
1248 if(bin_x1<1) bin_x1=1;
1249 Int_t bin_x2 = fXaxis.FindFixBin(x2);
1250 if(bin_x2>GetNbinsX()) bin_x2=GetNbinsX();
1251 Int_t bin_y1 = fYaxis.FindFixBin(y1);
1252 if(bin_y1<1) bin_y1=1;
1253 Int_t bin_y2 = fYaxis.FindFixBin(y2);
1254 if(bin_y2>GetNbinsY()) bin_y2=GetNbinsY();
1255 Int_t bin_q22 = GetBin(bin_x2,bin_y2);
1256 Int_t bin_q12 = GetBin(bin_x1,bin_y2);
1257 Int_t bin_q11 = GetBin(bin_x1,bin_y1);
1258 Int_t bin_q21 = GetBin(bin_x2,bin_y1);
1259 Double_t q11 = RetrieveBinContent(bin_q11);
1260 Double_t q12 = RetrieveBinContent(bin_q12);
1261 Double_t q21 = RetrieveBinContent(bin_q21);
1262 Double_t q22 = RetrieveBinContent(bin_q22);
1263 Double_t d = 1.0*(x2-x1)*(y2-y1);
1264 f = 1.0*q11/d*(x2-x)*(y2-y)+1.0*q21/d*(x-x1)*(y2-y)+1.0*q12/d*(x2-x)*(y-y1)+1.0*q22/d*(x-x1)*(y-y1);
1272 Double_t TH2::Interpolate(Double_t, Double_t, Double_t)
const
1274 Error(
"Interpolate",
"This function must be called with 2 arguments for a TH2");
1302 Double_t TH2::KolmogorovTest(
const TH1 *h2, Option_t *option)
const
1304 TString opt = option;
1308 TH1 *h1 = (TH1*)
this;
1309 if (h2 == 0)
return 0;
1310 const TAxis *xaxis1 = h1->GetXaxis();
1311 const TAxis *xaxis2 = h2->GetXaxis();
1312 const TAxis *yaxis1 = h1->GetYaxis();
1313 const TAxis *yaxis2 = h2->GetYaxis();
1314 Int_t ncx1 = xaxis1->GetNbins();
1315 Int_t ncx2 = xaxis2->GetNbins();
1316 Int_t ncy1 = yaxis1->GetNbins();
1317 Int_t ncy2 = yaxis2->GetNbins();
1320 if (h1->GetDimension() != 2 || h2->GetDimension() != 2) {
1321 Error(
"KolmogorovTest",
"Histograms must be 2-D\n");
1327 Error(
"KolmogorovTest",
"Number of channels in X is different, %d and %d\n",ncx1,ncx2);
1331 Error(
"KolmogorovTest",
"Number of channels in Y is different, %d and %d\n",ncy1,ncy2);
1336 Bool_t afunc1 = kFALSE;
1337 Bool_t afunc2 = kFALSE;
1338 Double_t difprec = 1e-5;
1339 Double_t diff1 = TMath::Abs(xaxis1->GetXmin() - xaxis2->GetXmin());
1340 Double_t diff2 = TMath::Abs(xaxis1->GetXmax() - xaxis2->GetXmax());
1341 if (diff1 > difprec || diff2 > difprec) {
1342 Error(
"KolmogorovTest",
"histograms with different binning along X");
1345 diff1 = TMath::Abs(yaxis1->GetXmin() - yaxis2->GetXmin());
1346 diff2 = TMath::Abs(yaxis1->GetXmax() - yaxis2->GetXmax());
1347 if (diff1 > difprec || diff2 > difprec) {
1348 Error(
"KolmogorovTest",
"histograms with different binning along Y");
1353 Int_t ibeg = 1, jbeg = 1;
1354 Int_t iend = ncx1, jend = ncy1;
1355 if (opt.Contains(
"U")) {ibeg = 0; jbeg = 0;}
1356 if (opt.Contains(
"O")) {iend = ncx1+1; jend = ncy1+1;}
1363 for (i = ibeg; i <= iend; i++) {
1364 for (j = jbeg; j <= jend; j++) {
1365 sum1 += h1->GetBinContent(i,j);
1366 sum2 += h2->GetBinContent(i,j);
1367 Double_t ew1 = h1->GetBinError(i,j);
1368 Double_t ew2 = h2->GetBinError(i,j);
1377 Error(
"KolmogorovTest",
"Integral is zero for h1=%s\n",h1->GetName());
1381 Error(
"KolmogorovTest",
"Integral is zero for h2=%s\n",h2->GetName());
1387 Double_t esum1 = 0, esum2 = 0;
1389 esum1 = sum1 * sum1 / w1;
1394 esum2 = sum2 * sum2 / w2;
1398 if (afunc2 && afunc1) {
1399 Error(
"KolmogorovTest",
"Errors are zero for both histograms\n");
1404 Double_t s1 = 1/sum1;
1405 Double_t s2 = 1/sum2;
1406 Double_t dfmax1 = 0;
1407 Double_t rsum1=0, rsum2=0;
1408 for (i=ibeg;i<=iend;i++) {
1409 for (j=jbeg;j<=jend;j++) {
1410 rsum1 += s1*h1->GetBinContent(i,j);
1411 rsum2 += s2*h2->GetBinContent(i,j);
1412 dfmax1 = TMath::Max(dfmax1, TMath::Abs(rsum1-rsum2));
1417 Double_t dfmax2 = 0;
1419 for (j=jbeg;j<=jend;j++) {
1420 for (i=ibeg;i<=iend;i++) {
1421 rsum1 += s1*h1->GetBinContent(i,j);
1422 rsum2 += s2*h2->GetBinContent(i,j);
1423 dfmax2 = TMath::Max(dfmax2, TMath::Abs(rsum1-rsum2));
1429 if (afunc1) factnm = TMath::Sqrt(esum2);
1430 else if (afunc2) factnm = TMath::Sqrt(esum1);
1431 else factnm = TMath::Sqrt(esum1*sum2/(esum1+esum2));
1434 Double_t dfmax = 0.5*(dfmax1+dfmax2);
1435 Double_t z = dfmax*factnm;
1437 prb = TMath::KolmogorovProb(z);
1439 Double_t prb1 = 0, prb2 = 0;
1441 if (opt.Contains(
"N") && !(afunc1 || afunc2 ) ) {
1444 Double_t d12 = esum1-esum2;
1445 Double_t chi2 = d12*d12/(esum1+esum2);
1446 prb2 = TMath::Prob(chi2,1);
1448 if (prb > 0 && prb2 > 0) prb = prb*prb2*(1-TMath::Log(prb*prb2));
1453 if (opt.Contains(
"D")) {
1454 printf(
" Kolmo Prob h1 = %s, sum1=%g\n",h1->GetName(),sum1);
1455 printf(
" Kolmo Prob h2 = %s, sum2=%g\n",h2->GetName(),sum2);
1456 printf(
" Kolmo Probabil = %f, Max Dist = %g\n",prb,dfmax);
1457 if (opt.Contains(
"N"))
1458 printf(
" Kolmo Probabil = %f for shape alone, =%f for normalisation alone\n",prb1,prb2);
1461 if (TMath::Abs(rsum1-1) > 0.002) Warning(
"KolmogorovTest",
"Numerical problems with h1=%s\n",h1->GetName());
1462 if (TMath::Abs(rsum2-1) > 0.002) Warning(
"KolmogorovTest",
"Numerical problems with h2=%s\n",h2->GetName());
1464 if(opt.Contains(
"M"))
return dfmax;
1474 TH2 *TH2::RebinX(Int_t ngroup,
const char *newname)
1476 return Rebin2D(ngroup, 1, newname);
1484 TH2 *TH2::RebinY(Int_t ngroup,
const char *newname)
1486 return Rebin2D(1, ngroup, newname);
1495 TH2 * TH2::Rebin( Int_t ngroup,
const char*newname,
const Double_t *xbins)
1497 if (xbins !=
nullptr) {
1498 Error(
"Rebin",
"Rebinning a 2-d histogram into variable bins is not supported (it is possible only for 1-d histograms). Return a nullptr");
1501 Info(
"Rebin",
"Rebinning only the x-axis. Use Rebin2D for rebinning both axes");
1502 return RebinX(ngroup, newname);
1530 TH2 *TH2::Rebin2D(Int_t nxgroup, Int_t nygroup,
const char *newname)
1532 Int_t nxbins = fXaxis.GetNbins();
1533 Int_t nybins = fYaxis.GetNbins();
1534 Int_t nx = nxbins + 2;
1535 Int_t ny = nybins + 2;
1536 Double_t xmin = fXaxis.GetXmin();
1537 Double_t xmax = fXaxis.GetXmax();
1538 Double_t ymin = fYaxis.GetXmin();
1539 Double_t ymax = fYaxis.GetXmax();
1541 if (GetDimension() != 2) {
1542 Error(
"Rebin2D",
"Histogram must be TH2. This histogram has %d dimensions.", GetDimension());
1545 if ((nxgroup <= 0) || (nxgroup > nxbins)) {
1546 Error(
"Rebin2D",
"Illegal value of nxgroup=%d",nxgroup);
1549 if ((nygroup <= 0) || (nygroup > nybins)) {
1550 Error(
"Rebin2D",
"Illegal value of nygroup=%d",nygroup);
1554 Int_t newxbins = nxbins / nxgroup;
1555 Int_t newybins = nybins / nygroup;
1556 Int_t newnx = newxbins + 2;
1557 Int_t newny = newybins + 2;
1560 Double_t *oldBins =
new Double_t[fNcells];
1561 for (Int_t i = 0; i < fNcells; ++i) oldBins[i] = RetrieveBinContent(i);
1563 Double_t* oldErrors = NULL;
1565 oldErrors =
new Double_t[fNcells];
1566 for (Int_t i = 0; i < fNcells; ++i) oldErrors[i] = GetBinErrorSqUnchecked(i);
1571 if (newname && strlen(newname)) {
1572 hnew = (TH2*)Clone();
1573 hnew->SetName(newname);
1576 bool resetStat =
false;
1579 if(newxbins * nxgroup != nxbins) {
1580 xmax = fXaxis.GetBinUpEdge(newxbins * nxgroup);
1583 if(newybins * nygroup != nybins) {
1584 ymax = fYaxis.GetBinUpEdge(newybins * nygroup);
1589 Int_t nXdivisions = fXaxis.GetNdivisions();
1590 Color_t xAxisColor = fXaxis.GetAxisColor();
1591 Color_t xLabelColor = fXaxis.GetLabelColor();
1592 Style_t xLabelFont = fXaxis.GetLabelFont();
1593 Float_t xLabelOffset = fXaxis.GetLabelOffset();
1594 Float_t xLabelSize = fXaxis.GetLabelSize();
1595 Float_t xTickLength = fXaxis.GetTickLength();
1596 Float_t xTitleOffset = fXaxis.GetTitleOffset();
1597 Float_t xTitleSize = fXaxis.GetTitleSize();
1598 Color_t xTitleColor = fXaxis.GetTitleColor();
1599 Style_t xTitleFont = fXaxis.GetTitleFont();
1601 Int_t nYdivisions = fYaxis.GetNdivisions();
1602 Color_t yAxisColor = fYaxis.GetAxisColor();
1603 Color_t yLabelColor = fYaxis.GetLabelColor();
1604 Style_t yLabelFont = fYaxis.GetLabelFont();
1605 Float_t yLabelOffset = fYaxis.GetLabelOffset();
1606 Float_t yLabelSize = fYaxis.GetLabelSize();
1607 Float_t yTickLength = fYaxis.GetTickLength();
1608 Float_t yTitleOffset = fYaxis.GetTitleOffset();
1609 Float_t yTitleSize = fYaxis.GetTitleSize();
1610 Color_t yTitleColor = fYaxis.GetTitleColor();
1611 Style_t yTitleFont = fYaxis.GetTitleFont();
1615 if (nxgroup != 1 || nygroup != 1) {
1616 if(fXaxis.GetXbins()->GetSize() > 0 || fYaxis.GetXbins()->GetSize() > 0){
1618 Double_t *xbins =
new Double_t[newxbins + 1];
1619 for(Int_t i = 0; i <= newxbins; ++i) xbins[i] = fXaxis.GetBinLowEdge(1 + i * nxgroup);
1620 Double_t *ybins =
new Double_t[newybins + 1];
1621 for(Int_t i = 0; i <= newybins; ++i) ybins[i] = fYaxis.GetBinLowEdge(1 + i * nygroup);
1622 hnew->SetBins(newxbins, xbins, newybins, ybins);
1626 hnew->SetBins(newxbins, xmin, xmax, newybins, ymin, ymax);
1630 hnew->UpdateBinContent(0, oldBins[0]);
1631 if (oldErrors) hnew->fSumw2[0] = 0;
1634 for(Int_t binx = 1, oldbinx = 1; binx < newnx; ++binx, oldbinx += nxgroup){
1635 Double_t binContent = 0.0, binErrorSq = 0.0;
1636 for (Int_t i = 0; i < nxgroup && (oldbinx + i) < nx; ++i) {
1637 Int_t bin = oldbinx + i;
1638 binContent += oldBins[bin];
1639 if(oldErrors) binErrorSq += oldErrors[bin];
1641 Int_t newbin = binx;
1642 hnew->UpdateBinContent(newbin, binContent);
1643 if (oldErrors) hnew->fSumw2[newbin] = binErrorSq;
1647 for(Int_t biny = 1, oldbiny = 1; biny < newny; ++biny, oldbiny += nygroup){
1648 Double_t binContent = 0.0, binErrorSq = 0.0;
1649 for (Int_t j = 0; j < nygroup && (oldbiny + j) < ny; ++j) {
1650 Int_t bin = (oldbiny + j) * nx;
1651 binContent += oldBins[bin];
1652 if(oldErrors) binErrorSq += oldErrors[bin];
1654 Int_t newbin = biny * newnx;
1655 hnew->UpdateBinContent(newbin, binContent);
1656 if (oldErrors) hnew->fSumw2[newbin] = binErrorSq;
1660 for (Int_t binx = 1, oldbinx = 1; binx < newnx; ++binx, oldbinx += nxgroup) {
1661 for (Int_t biny = 1, oldbiny = 1; biny < newny; ++biny, oldbiny += nygroup) {
1662 Double_t binContent = 0.0, binErrorSq = 0.0;
1663 for (Int_t i = 0; i < nxgroup && (oldbinx + i) < nx; ++i) {
1664 for (Int_t j = 0; j < nygroup && (oldbiny + j) < ny; ++j) {
1665 Int_t bin = oldbinx + i + (oldbiny + j) * nx;
1666 binContent += oldBins[bin];
1667 if (oldErrors) binErrorSq += oldErrors[bin];
1670 Int_t newbin = binx + biny * newnx;
1671 hnew->UpdateBinContent(newbin, binContent);
1672 if (oldErrors) hnew->fSumw2[newbin] = binErrorSq;
1678 fXaxis.SetNdivisions(nXdivisions);
1679 fXaxis.SetAxisColor(xAxisColor);
1680 fXaxis.SetLabelColor(xLabelColor);
1681 fXaxis.SetLabelFont(xLabelFont);
1682 fXaxis.SetLabelOffset(xLabelOffset);
1683 fXaxis.SetLabelSize(xLabelSize);
1684 fXaxis.SetTickLength(xTickLength);
1685 fXaxis.SetTitleOffset(xTitleOffset);
1686 fXaxis.SetTitleSize(xTitleSize);
1687 fXaxis.SetTitleColor(xTitleColor);
1688 fXaxis.SetTitleFont(xTitleFont);
1690 fYaxis.SetNdivisions(nYdivisions);
1691 fYaxis.SetAxisColor(yAxisColor);
1692 fYaxis.SetLabelColor(yLabelColor);
1693 fYaxis.SetLabelFont(yLabelFont);
1694 fYaxis.SetLabelOffset(yLabelOffset);
1695 fYaxis.SetLabelSize(yLabelSize);
1696 fYaxis.SetTickLength(yTickLength);
1697 fYaxis.SetTitleOffset(yTitleOffset);
1698 fYaxis.SetTitleSize(yTitleSize);
1699 fYaxis.SetTitleColor(yTitleColor);
1700 fYaxis.SetTitleFont(yTitleFont);
1702 if (resetStat) hnew->ResetStats();
1705 if (oldErrors)
delete [] oldErrors;
1712 TProfile *TH2::DoProfile(
bool onX,
const char *name, Int_t firstbin, Int_t lastbin, Option_t *option)
const
1714 TString opt = option;
1717 Int_t i1 = opt.Index(
"[");
1719 Int_t i2 = opt.Index(
"]");
1720 cut = opt(i1,i2-i1+1);
1723 bool originalRange = opt.Contains(
"o");
1725 const TAxis& outAxis = ( onX ? fXaxis : fYaxis );
1726 const TAxis& inAxis = ( onX ? fYaxis : fXaxis );
1727 Int_t inN = inAxis.GetNbins();
1728 const char *expectedName = ( onX ?
"_pfx" :
"_pfy" );
1732 Int_t firstOutBin = std::max(outAxis.GetFirst(),1);
1733 Int_t lastOutBin = std::min(outAxis.GetLast(),outAxis.GetNbins() ) ;
1735 if ( lastbin < firstbin && inAxis.TestBit(TAxis::kAxisRange) ) {
1736 firstbin = inAxis.GetFirst();
1737 lastbin = inAxis.GetLast();
1741 if (firstbin == 0 && lastbin == 0)
1744 lastbin = inAxis.GetNbins();
1747 if (firstbin < 0) firstbin = 1;
1748 if (lastbin < 0) lastbin = inN;
1749 if (lastbin > inN+1) lastbin = inN;
1752 char *pname = (
char*)name;
1753 if (name && strcmp(name, expectedName) == 0) {
1754 Int_t nch = strlen(GetName()) + 5;
1755 pname =
new char[nch];
1756 snprintf(pname,nch,
"%s%s",GetName(),name);
1761 TObject *h1obj = gROOT->FindObject(pname);
1762 if (h1obj && h1obj->InheritsFrom(TH1::Class())) {
1763 if (h1obj->IsA() != TProfile::Class() ) {
1764 Error(
"DoProfile",
"Histogram with name %s must be a TProfile and is a %s",name,h1obj->ClassName());
1767 h1 = (TProfile*)h1obj;
1772 const TArrayD *xbins = outAxis.GetXbins();
1773 if (xbins->fN == 0) {
1774 if ( originalRange )
1775 h1->SetBins(outAxis.GetNbins(),outAxis.GetXmin(),outAxis.GetXmax());
1777 h1->SetBins(lastOutBin-firstOutBin+1,outAxis.GetBinLowEdge(firstOutBin),outAxis.GetBinUpEdge(lastOutBin));
1781 h1->SetBins(outAxis.GetNbins(),xbins->fArray);
1783 h1->SetBins(lastOutBin-firstOutBin+1,&xbins->fArray[firstOutBin-1]);
1788 if (opt.Contains(
"[")) {
1789 ((TH2 *)
this)->GetPainter();
1790 if (fPainter) ncuts = fPainter->MakeCuts((
char*)cut.Data());
1794 const TArrayD *bins = outAxis.GetXbins();
1795 if (bins->fN == 0) {
1796 if ( originalRange )
1797 h1 =
new TProfile(pname,GetTitle(),outAxis.GetNbins(),outAxis.GetXmin(),outAxis.GetXmax(),opt);
1799 h1 =
new TProfile(pname,GetTitle(),lastOutBin-firstOutBin+1,
1800 outAxis.GetBinLowEdge(firstOutBin),
1801 outAxis.GetBinUpEdge(lastOutBin), opt);
1805 h1 =
new TProfile(pname,GetTitle(),outAxis.GetNbins(),bins->fArray,opt);
1807 h1 =
new TProfile(pname,GetTitle(),lastOutBin-firstOutBin+1,&bins->fArray[firstOutBin-1],opt);
1810 if (pname != name)
delete [] pname;
1813 h1->GetXaxis()->ImportAttributes( &outAxis);
1814 h1->SetLineColor(this->GetLineColor());
1815 h1->SetFillColor(this->GetFillColor());
1816 h1->SetMarkerColor(this->GetMarkerColor());
1817 h1->SetMarkerStyle(this->GetMarkerStyle());
1821 TArrayD & binSumw2 = *(h1->GetBinSumw2());
1822 bool useWeights = (GetSumw2N() > 0);
1823 if (useWeights && (binSumw2.fN != h1->GetNcells()) ) h1->Sumw2();
1826 else h1->SetBit(TH1::kIsNotW);
1830 Double_t totcont = 0;
1835 for ( Int_t outbin = 0; outbin <= outAxis.GetNbins() + 1; ++outbin) {
1836 if (outAxis.TestBit(TAxis::kAxisRange) && ( outbin < firstOutBin || outbin > lastOutBin ))
continue;
1839 Double_t xOut = outAxis.GetBinCenter(outbin);
1840 Int_t binOut = h1->GetXaxis()->FindBin( xOut );
1841 if (binOut <0)
continue;
1843 for (Int_t inbin = firstbin ; inbin <= lastbin ; ++inbin) {
1845 if (onX) { binx = outbin; biny=inbin; }
1846 else { binx = inbin; biny=outbin; }
1849 if (!fPainter->IsInside(binx,biny))
continue;
1851 Int_t bin = GetBin(binx, biny);
1852 Double_t cxy = RetrieveBinContent(bin);
1858 if ( useWeights ) tmp = binSumw2.fArray[binOut];
1859 h1->Fill( xOut, inAxis.GetBinCenter(inbin), cxy );
1860 if ( useWeights ) binSumw2.fArray[binOut] = tmp + fSumw2.fArray[bin];
1873 h1->SetEntries( h1->GetEffectiveEntries() );
1876 if (opt.Contains(
"d")) {
1877 TVirtualPad *padsav = gPad;
1878 TVirtualPad *pad = gROOT->GetSelectedPad();
1880 opt.Remove(opt.First(
"d"),1);
1881 if (!gPad || !gPad->FindObject(h1)) {
1886 if (padsav) padsav->cd();
1935 TProfile *TH2::ProfileX(
const char *name, Int_t firstybin, Int_t lastybin, Option_t *option)
const
1937 return DoProfile(
true, name, firstybin, lastybin, option);
1985 TProfile *TH2::ProfileY(
const char *name, Int_t firstxbin, Int_t lastxbin, Option_t *option)
const
1987 return DoProfile(
false, name, firstxbin, lastxbin, option);
1995 TH1D *TH2::DoProjection(
bool onX,
const char *name, Int_t firstbin, Int_t lastbin, Option_t *option)
const
1997 const char *expectedName = 0;
1999 const TAxis* outAxis;
2000 const TAxis* inAxis;
2002 TString opt = option;
2004 Int_t i1 = opt.Index(
"[");
2006 Int_t i2 = opt.Index(
"]");
2007 cut = opt(i1,i2-i1+1);
2010 bool originalRange = opt.Contains(
"o");
2014 expectedName =
"_px";
2015 inNbin = fYaxis.GetNbins();
2016 outAxis = GetXaxis();
2017 inAxis = GetYaxis();
2021 expectedName =
"_py";
2022 inNbin = fXaxis.GetNbins();
2023 outAxis = GetYaxis();
2024 inAxis = GetXaxis();
2029 Int_t firstOutBin = std::max(outAxis->GetFirst(),1);
2030 Int_t lastOutBin = std::min(outAxis->GetLast(),outAxis->GetNbins() ) ;
2032 if ( lastbin < firstbin && inAxis->TestBit(TAxis::kAxisRange) ) {
2033 firstbin = inAxis->GetFirst();
2034 lastbin = inAxis->GetLast();
2038 if (firstbin == 0 && lastbin == 0)
2041 lastbin = inAxis->GetNbins();
2044 if (firstbin < 0) firstbin = 0;
2045 if (lastbin < 0) lastbin = inNbin + 1;
2046 if (lastbin > inNbin+1) lastbin = inNbin + 1;
2049 char *pname = (
char*)name;
2050 if (name && strcmp(name,expectedName) == 0) {
2051 Int_t nch = strlen(GetName()) + 4;
2052 pname =
new char[nch];
2053 snprintf(pname,nch,
"%s%s",GetName(),name);
2059 TObject *h1obj = gROOT->FindObject(pname);
2060 if (h1obj && h1obj->InheritsFrom(TH1::Class())) {
2061 if (h1obj->IsA() != TH1D::Class() ) {
2062 Error(
"DoProjection",
"Histogram with name %s must be a TH1D and is a %s",name,h1obj->ClassName());
2070 const TArrayD *xbins = outAxis->GetXbins();
2071 if (xbins->fN == 0) {
2072 if ( originalRange )
2073 h1->SetBins(outAxis->GetNbins(),outAxis->GetXmin(),outAxis->GetXmax());
2075 h1->SetBins(lastOutBin-firstOutBin+1,outAxis->GetBinLowEdge(firstOutBin),outAxis->GetBinUpEdge(lastOutBin));
2079 h1->SetBins(outAxis->GetNbins(),xbins->fArray);
2081 h1->SetBins(lastOutBin-firstOutBin+1,&xbins->fArray[firstOutBin-1]);
2086 if (opt.Contains(
"[")) {
2087 ((TH2 *)
this)->GetPainter();
2088 if (fPainter) ncuts = fPainter->MakeCuts((
char*)cut.Data());
2092 const TArrayD *bins = outAxis->GetXbins();
2093 if (bins->fN == 0) {
2094 if ( originalRange )
2095 h1 =
new TH1D(pname,GetTitle(),outAxis->GetNbins(),outAxis->GetXmin(),outAxis->GetXmax());
2097 h1 =
new TH1D(pname,GetTitle(),lastOutBin-firstOutBin+1,
2098 outAxis->GetBinLowEdge(firstOutBin),outAxis->GetBinUpEdge(lastOutBin));
2102 h1 =
new TH1D(pname,GetTitle(),outAxis->GetNbins(),bins->fArray);
2104 h1 =
new TH1D(pname,GetTitle(),lastOutBin-firstOutBin+1,&bins->fArray[firstOutBin-1]);
2106 if (opt.Contains(
"e") || GetSumw2N() ) h1->Sumw2();
2108 if (pname != name)
delete [] pname;
2111 h1->GetXaxis()->ImportAttributes(outAxis);
2112 THashList* labels=outAxis->GetLabels();
2117 while ((lb=(TObjString*)iL())) {
2118 h1->GetXaxis()->SetBinLabel(i,lb->String().Data());
2123 h1->SetLineColor(this->GetLineColor());
2124 h1->SetFillColor(this->GetFillColor());
2125 h1->SetMarkerColor(this->GetMarkerColor());
2126 h1->SetMarkerStyle(this->GetMarkerStyle());
2130 Double_t totcont = 0;
2131 Bool_t computeErrors = h1->GetSumw2N();
2136 for ( Int_t outbin = 0; outbin <= outAxis->GetNbins() + 1; ++outbin) {
2139 if (outAxis->TestBit(TAxis::kAxisRange) && ( outbin < firstOutBin || outbin > lastOutBin ))
continue;
2141 for (Int_t inbin = firstbin ; inbin <= lastbin ; ++inbin) {
2143 if (onX) { binx = outbin; biny=inbin; }
2144 else { binx = inbin; biny=outbin; }
2147 if (!fPainter->IsInside(binx,biny))
continue;
2150 cont += GetBinContent(binx,biny);
2151 if (computeErrors) {
2152 Double_t exy = GetBinError(binx,biny);
2157 Int_t binOut = h1->GetXaxis()->FindBin( outAxis->GetBinCenter(outbin) );
2158 h1->SetBinContent(binOut ,cont);
2159 if (computeErrors) h1->SetBinError(binOut,TMath::Sqrt(err2));
2165 bool reuseStats =
false;
2166 if ( ( GetStatOverflowsBehaviour() ==
false && firstbin == 1 && lastbin == inNbin ) ||
2167 ( GetStatOverflowsBehaviour() ==
true && firstbin == 0 && lastbin == inNbin + 1 ) )
2171 double eps = 1.E-12;
2172 if (IsA() == TH2F::Class() ) eps = 1.E-6;
2173 if (fTsumw != 0 && TMath::Abs( fTsumw - totcont) < TMath::Abs(fTsumw) * eps)
2176 if (ncuts) reuseStats =
false;
2178 bool reuseEntries = reuseStats;
2180 reuseEntries &= (firstbin==0 && lastbin == inNbin+1);
2182 Double_t stats[kNstat];
2185 stats[2] = stats[4];
2186 stats[3] = stats[5];
2188 h1->PutStats(stats);
2194 h1->SetEntries( h1->GetEffectiveEntries() );
2197 h1->SetEntries(fEntries);
2204 Double_t entries = TMath::Floor( totcont + 0.5);
2205 if (h1->GetSumw2N()) entries = h1->GetEffectiveEntries();
2206 h1->SetEntries( entries );
2209 if (opt.Contains(
"d")) {
2210 TVirtualPad *padsav = gPad;
2211 TVirtualPad *pad = gROOT->GetSelectedPad();
2213 opt.Remove(opt.First(
"d"),1);
2215 if (opt.Contains(
"e")) opt.Remove(opt.First(
"e"),1);
2216 if (!gPad || !gPad->FindObject(h1)) {
2221 if (padsav) padsav->cd();
2261 TH1D *TH2::ProjectionX(
const char *name, Int_t firstybin, Int_t lastybin, Option_t *option)
const
2263 return DoProjection(
true, name, firstybin, lastybin, option);
2301 TH1D *TH2::ProjectionY(
const char *name, Int_t firstxbin, Int_t lastxbin, Option_t *option)
const
2303 return DoProjection(
false, name, firstxbin, lastxbin, option);
2310 void TH2::PutStats(Double_t *stats)
2312 TH1::PutStats(stats);
2314 fTsumwy2 = stats[5];
2315 fTsumwxy = stats[6];
2330 TH1D* TH2::QuantilesX( Double_t prob,
const char * name)
const
2332 return DoQuantiles(
true, name, prob);
2343 TH1D* TH2::QuantilesY( Double_t prob,
const char * name)
const
2345 return DoQuantiles(
false, name, prob);
2352 TH1D* TH2::DoQuantiles(
bool onX,
const char * name, Double_t prob)
const
2354 const TAxis *outAxis = 0;
2356 outAxis = GetXaxis();
2358 outAxis = GetYaxis();
2362 TString qname = name;
2363 if (qname.IsNull() || qname ==
"_qx" || qname ==
"_qy") {
2364 const char * qtype = (onX) ?
"qx" :
"qy";
2365 qname = TString::Format(
"%s_%s_%3.2f",GetName(),qtype, prob);
2370 TObject *h1obj = gROOT->FindObject(qname);
2372 h1 =
dynamic_cast<TH1D*
>(h1obj);
2374 Error(
"DoQuantiles",
"Histogram with name %s must be a TH1D and is a %s",qname.Data(),h1obj->ClassName());
2382 h1 =
new TH1D(qname, GetTitle(), 1, 0, 1);
2385 Int_t firstOutBin = std::max(outAxis->GetFirst(),1);
2386 Int_t lastOutBin = std::max(outAxis->GetLast(),outAxis->GetNbins());
2387 const TArrayD *xbins = outAxis->GetXbins();
2389 h1->SetBins(lastOutBin-firstOutBin+1,outAxis->GetBinLowEdge(firstOutBin),outAxis->GetBinUpEdge(lastOutBin));
2391 h1->SetBins(lastOutBin-firstOutBin+1,&xbins->fArray[firstOutBin-1]);
2398 for (
int ibin = outAxis->GetFirst() ; ibin <= outAxis->GetLast() ; ++ibin) {
2401 slice = DoProjection(!onX,
"tmp",ibin,ibin,
"");
2403 if (slice->GetSum() == 0)
continue;
2404 slice->GetQuantiles(1,qq,pp);
2405 h1->SetBinContent(ibin,qq[0]);
2409 Double_t n = slice->GetEffectiveEntries();
2410 Double_t f = TMath::Gaus(qq[0], slice->GetMean(), slice->GetStdDev(), kTRUE);
2414 error = TMath::Sqrt( prob*(1.-prob)/ (n * f * f) );
2415 h1->SetBinError(ibin, error);
2417 if (slice)
delete slice;
2425 void TH2::Reset(Option_t *option)
2428 TString opt = option;
2431 if (opt.Contains(
"ICE") && !opt.Contains(
"S"))
return;
2441 void TH2::SetBinContent(Int_t bin, Double_t content)
2445 if (bin < 0)
return;
2446 if (bin >= fNcells)
return;
2447 UpdateBinContent(bin, content);
2458 void TH2::SetShowProjectionX(Int_t nbins)
2462 if (fPainter) fPainter->SetShowProjection(
"x",nbins);
2473 void TH2::SetShowProjectionY(Int_t nbins)
2477 if (fPainter) fPainter->SetShowProjection(
"y",nbins);
2486 TH1 *TH2::ShowBackground(Int_t niter, Option_t *option)
2489 return (TH1*)gROOT->ProcessLineFast(Form(
"TSpectrum2::StaticBackground((TH1*)0x%lx,%d,\"%s\")",
2490 (ULong_t)
this, niter, option));
2502 Int_t TH2::ShowPeaks(Double_t sigma, Option_t *option, Double_t threshold)
2505 return (Int_t)gROOT->ProcessLineFast(Form(
"TSpectrum2::StaticSearch((TH1*)0x%lx,%g,\"%s\",%g)",
2506 (ULong_t)
this, sigma, option, threshold));
2531 void TH2::Smooth(Int_t ntimes, Option_t *option)
2533 Double_t k5a[5][5] = { { 0, 0, 1, 0, 0 },
2537 { 0, 0, 1, 0, 0 } };
2538 Double_t k5b[5][5] = { { 0, 1, 2, 1, 0 },
2542 { 0, 1, 2, 1, 0 } };
2543 Double_t k3a[3][3] = { { 0, 1, 0 },
2548 Warning(
"Smooth",
"Currently only ntimes=1 is supported");
2550 TString opt = option;
2554 Double_t *kernel = &k5a[0][0];
2555 if (opt.Contains(
"k5b")) kernel = &k5b[0][0];
2556 if (opt.Contains(
"k3a")) {
2557 kernel = &k3a[0][0];
2563 Int_t ifirst = fXaxis.GetFirst();
2564 Int_t ilast = fXaxis.GetLast();
2565 Int_t jfirst = fYaxis.GetFirst();
2566 Int_t jlast = fYaxis.GetLast();
2569 Double_t nentries = fEntries;
2570 Int_t nx = GetNbinsX();
2571 Int_t ny = GetNbinsY();
2572 Int_t bufSize = (nx+2)*(ny+2);
2573 Double_t *buf =
new Double_t[bufSize];
2575 if (fSumw2.fN) ebuf =
new Double_t[bufSize];
2579 for (i=ifirst; i<=ilast; i++){
2580 for (j=jfirst; j<=jlast; j++){
2582 buf[bin] = RetrieveBinContent(bin);
2583 if (ebuf) ebuf[bin]=GetBinError(bin);
2588 Int_t x_push = (ksize_x-1)/2;
2589 Int_t y_push = (ksize_y-1)/2;
2592 for (i=ifirst; i<=ilast; i++){
2593 for (j=jfirst; j<=jlast; j++) {
2594 Double_t content = 0.0;
2595 Double_t error = 0.0;
2596 Double_t norm = 0.0;
2598 for (Int_t n=0; n<ksize_x; n++) {
2599 for (Int_t m=0; m<ksize_y; m++) {
2600 Int_t xb = i+(n-x_push);
2601 Int_t yb = j+(m-y_push);
2602 if ( (xb >= 1) && (xb <= nx) && (yb >= 1) && (yb <= ny) ) {
2603 bin = GetBin(xb,yb);
2604 Double_t k = kernel[n*ksize_y +m];
2608 content += k*buf[bin];
2609 if (ebuf) error += k*k*ebuf[bin]*ebuf[bin];
2615 if ( norm != 0.0 ) {
2616 SetBinContent(i,j,content/norm);
2618 error /= (norm*norm);
2619 SetBinError(i,j,sqrt(error));
2624 fEntries = nentries;
2634 void TH2::Streamer(TBuffer &R__b)
2636 if (R__b.IsReading()) {
2638 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
2640 R__b.ReadClassBuffer(TH2::Class(),
this, R__v, R__s, R__c);
2644 TH1::Streamer(R__b);
2645 R__b >> fScalefactor;
2652 R__b.WriteClassBuffer(TH2::Class(),
this);
2668 TH2C::TH2C(): TH2(), TArrayC()
2671 if (fgDefaultSumw2) Sumw2();
2686 TH2C::TH2C(
const char *name,
const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
2687 ,Int_t nbinsy,Double_t ylow,Double_t yup)
2688 :TH2(name,title,nbinsx,xlow,xup,nbinsy,ylow,yup)
2690 TArrayC::Set(fNcells);
2691 if (fgDefaultSumw2) Sumw2();
2693 if (xlow >= xup || ylow >= yup) SetBuffer(fgBufferSize);
2700 TH2C::TH2C(
const char *name,
const char *title,Int_t nbinsx,
const Double_t *xbins
2701 ,Int_t nbinsy,Double_t ylow,Double_t yup)
2702 :TH2(name,title,nbinsx,xbins,nbinsy,ylow,yup)
2704 TArrayC::Set(fNcells);
2705 if (fgDefaultSumw2) Sumw2();
2712 TH2C::TH2C(
const char *name,
const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
2713 ,Int_t nbinsy,
const Double_t *ybins)
2714 :TH2(name,title,nbinsx,xlow,xup,nbinsy,ybins)
2716 TArrayC::Set(fNcells);
2717 if (fgDefaultSumw2) Sumw2();
2724 TH2C::TH2C(
const char *name,
const char *title,Int_t nbinsx,
const Double_t *xbins
2725 ,Int_t nbinsy,
const Double_t *ybins)
2726 :TH2(name,title,nbinsx,xbins,nbinsy,ybins)
2728 TArrayC::Set(fNcells);
2729 if (fgDefaultSumw2) Sumw2();
2736 TH2C::TH2C(
const char *name,
const char *title,Int_t nbinsx,
const Float_t *xbins
2737 ,Int_t nbinsy,
const Float_t *ybins)
2738 :TH2(name,title,nbinsx,xbins,nbinsy,ybins)
2740 TArrayC::Set(fNcells);
2741 if (fgDefaultSumw2) Sumw2();
2748 TH2C::TH2C(
const TH2C &h2c) : TH2(), TArrayC()
2750 ((TH2C&)h2c).Copy(*
this);
2757 void TH2C::AddBinContent(Int_t bin)
2759 if (fArray[bin] < 127) fArray[bin]++;
2766 void TH2C::AddBinContent(Int_t bin, Double_t w)
2768 Int_t newval = fArray[bin] + Int_t(w);
2769 if (newval > -128 && newval < 128) {fArray[bin] = Char_t(newval);
return;}
2770 if (newval < -127) fArray[bin] = -127;
2771 if (newval > 127) fArray[bin] = 127;
2778 void TH2C::Copy(TObject &newth2)
const
2780 TH2::Copy((TH2C&)newth2);
2787 void TH2C::Reset(Option_t *option)
2798 void TH2C::SetBinsLength(Int_t n)
2800 if (n < 0) n = (fXaxis.GetNbins()+2)*(fYaxis.GetNbins()+2);
2809 void TH2C::Streamer(TBuffer &R__b)
2811 if (R__b.IsReading()) {
2813 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
2815 R__b.ReadClassBuffer(TH2C::Class(),
this, R__v, R__s, R__c);
2821 TH1::Streamer(R__b);
2822 TArrayC::Streamer(R__b);
2824 R__b >> fScalefactor;
2829 TH2::Streamer(R__b);
2830 TArrayC::Streamer(R__b);
2831 R__b.CheckByteCount(R__s, R__c, TH2C::IsA());
2836 R__b.WriteClassBuffer(TH2C::Class(),
this);
2844 TH2C& TH2C::operator=(
const TH2C &h1)
2846 if (
this != &h1) ((TH2C&)h1).Copy(*
this);
2854 TH2C operator*(Float_t c1, TH2C &h1)
2858 hnew.SetDirectory(0);
2866 TH2C operator+(TH2C &h1, TH2C &h2)
2870 hnew.SetDirectory(0);
2878 TH2C operator-(TH2C &h1, TH2C &h2)
2882 hnew.SetDirectory(0);
2890 TH2C operator*(TH2C &h1, TH2C &h2)
2894 hnew.SetDirectory(0);
2902 TH2C operator/(TH2C &h1, TH2C &h2)
2906 hnew.SetDirectory(0);
2922 TH2S::TH2S(): TH2(), TArrayS()
2925 if (fgDefaultSumw2) Sumw2();
2940 TH2S::TH2S(
const char *name,
const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
2941 ,Int_t nbinsy,Double_t ylow,Double_t yup)
2942 :TH2(name,title,nbinsx,xlow,xup,nbinsy,ylow,yup)
2944 TArrayS::Set(fNcells);
2945 if (fgDefaultSumw2) Sumw2();
2947 if (xlow >= xup || ylow >= yup) SetBuffer(fgBufferSize);
2954 TH2S::TH2S(
const char *name,
const char *title,Int_t nbinsx,
const Double_t *xbins
2955 ,Int_t nbinsy,Double_t ylow,Double_t yup)
2956 :TH2(name,title,nbinsx,xbins,nbinsy,ylow,yup)
2958 TArrayS::Set(fNcells);
2959 if (fgDefaultSumw2) Sumw2();
2966 TH2S::TH2S(
const char *name,
const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
2967 ,Int_t nbinsy,
const Double_t *ybins)
2968 :TH2(name,title,nbinsx,xlow,xup,nbinsy,ybins)
2970 TArrayS::Set(fNcells);
2971 if (fgDefaultSumw2) Sumw2();
2978 TH2S::TH2S(
const char *name,
const char *title,Int_t nbinsx,
const Double_t *xbins
2979 ,Int_t nbinsy,
const Double_t *ybins)
2980 :TH2(name,title,nbinsx,xbins,nbinsy,ybins)
2982 TArrayS::Set(fNcells);
2983 if (fgDefaultSumw2) Sumw2();
2990 TH2S::TH2S(
const char *name,
const char *title,Int_t nbinsx,
const Float_t *xbins
2991 ,Int_t nbinsy,
const Float_t *ybins)
2992 :TH2(name,title,nbinsx,xbins,nbinsy,ybins)
2994 TArrayS::Set(fNcells);
2995 if (fgDefaultSumw2) Sumw2();
3002 TH2S::TH2S(
const TH2S &h2s) : TH2(), TArrayS()
3004 ((TH2S&)h2s).Copy(*
this);
3011 void TH2S::AddBinContent(Int_t bin)
3013 if (fArray[bin] < 32767) fArray[bin]++;
3020 void TH2S::AddBinContent(Int_t bin, Double_t w)
3022 Int_t newval = fArray[bin] + Int_t(w);
3023 if (newval > -32768 && newval < 32768) {fArray[bin] = Short_t(newval);
return;}
3024 if (newval < -32767) fArray[bin] = -32767;
3025 if (newval > 32767) fArray[bin] = 32767;
3032 void TH2S::Copy(TObject &newth2)
const
3034 TH2::Copy((TH2S&)newth2);
3041 void TH2S::Reset(Option_t *option)
3052 void TH2S::SetBinsLength(Int_t n)
3054 if (n < 0) n = (fXaxis.GetNbins()+2)*(fYaxis.GetNbins()+2);
3063 void TH2S::Streamer(TBuffer &R__b)
3065 if (R__b.IsReading()) {
3067 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
3069 R__b.ReadClassBuffer(TH2S::Class(),
this, R__v, R__s, R__c);
3075 TH1::Streamer(R__b);
3076 TArrayS::Streamer(R__b);
3078 R__b >> fScalefactor;
3083 TH2::Streamer(R__b);
3084 TArrayS::Streamer(R__b);
3085 R__b.CheckByteCount(R__s, R__c, TH2S::IsA());
3090 R__b.WriteClassBuffer(TH2S::Class(),
this);
3098 TH2S& TH2S::operator=(
const TH2S &h1)
3100 if (
this != &h1) ((TH2S&)h1).Copy(*
this);
3108 TH2S operator*(Float_t c1, TH2S &h1)
3112 hnew.SetDirectory(0);
3120 TH2S operator+(TH2S &h1, TH2S &h2)
3124 hnew.SetDirectory(0);
3132 TH2S operator-(TH2S &h1, TH2S &h2)
3136 hnew.SetDirectory(0);
3144 TH2S operator*(TH2S &h1, TH2S &h2)
3148 hnew.SetDirectory(0);
3156 TH2S operator/(TH2S &h1, TH2S &h2)
3160 hnew.SetDirectory(0);
3176 TH2I::TH2I(): TH2(), TArrayI()
3179 if (fgDefaultSumw2) Sumw2();
3194 TH2I::TH2I(
const char *name,
const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
3195 ,Int_t nbinsy,Double_t ylow,Double_t yup)
3196 :TH2(name,title,nbinsx,xlow,xup,nbinsy,ylow,yup)
3198 TArrayI::Set(fNcells);
3199 if (fgDefaultSumw2) Sumw2();
3201 if (xlow >= xup || ylow >= yup) SetBuffer(fgBufferSize);
3208 TH2I::TH2I(
const char *name,
const char *title,Int_t nbinsx,
const Double_t *xbins
3209 ,Int_t nbinsy,Double_t ylow,Double_t yup)
3210 :TH2(name,title,nbinsx,xbins,nbinsy,ylow,yup)
3212 TArrayI::Set(fNcells);
3213 if (fgDefaultSumw2) Sumw2();
3220 TH2I::TH2I(
const char *name,
const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
3221 ,Int_t nbinsy,
const Double_t *ybins)
3222 :TH2(name,title,nbinsx,xlow,xup,nbinsy,ybins)
3224 TArrayI::Set(fNcells);
3225 if (fgDefaultSumw2) Sumw2();
3232 TH2I::TH2I(
const char *name,
const char *title,Int_t nbinsx,
const Double_t *xbins
3233 ,Int_t nbinsy,
const Double_t *ybins)
3234 :TH2(name,title,nbinsx,xbins,nbinsy,ybins)
3236 TArrayI::Set(fNcells);
3237 if (fgDefaultSumw2) Sumw2();
3244 TH2I::TH2I(
const char *name,
const char *title,Int_t nbinsx,
const Float_t *xbins
3245 ,Int_t nbinsy,
const Float_t *ybins)
3246 :TH2(name,title,nbinsx,xbins,nbinsy,ybins)
3248 TArrayI::Set(fNcells);
3249 if (fgDefaultSumw2) Sumw2();
3256 TH2I::TH2I(
const TH2I &h2i) : TH2(), TArrayI()
3258 ((TH2I&)h2i).Copy(*
this);
3265 void TH2I::AddBinContent(Int_t bin)
3267 if (fArray[bin] < 2147483647) fArray[bin]++;
3274 void TH2I::AddBinContent(Int_t bin, Double_t w)
3276 Long64_t newval = fArray[bin] + Long64_t(w);
3277 if (newval > -2147483647 && newval < 2147483647) {fArray[bin] = Int_t(newval);
return;}
3278 if (newval < -2147483647) fArray[bin] = -2147483647;
3279 if (newval > 2147483647) fArray[bin] = 2147483647;
3286 void TH2I::Copy(TObject &newth2)
const
3288 TH2::Copy((TH2I&)newth2);
3295 void TH2I::Reset(Option_t *option)
3306 void TH2I::SetBinsLength(Int_t n)
3308 if (n < 0) n = (fXaxis.GetNbins()+2)*(fYaxis.GetNbins()+2);
3317 TH2I& TH2I::operator=(
const TH2I &h1)
3319 if (
this != &h1) ((TH2I&)h1).Copy(*
this);
3327 TH2I operator*(Float_t c1, TH2I &h1)
3331 hnew.SetDirectory(0);
3339 TH2I operator+(TH2I &h1, TH2I &h2)
3343 hnew.SetDirectory(0);
3351 TH2I operator-(TH2I &h1, TH2I &h2)
3355 hnew.SetDirectory(0);
3363 TH2I operator*(TH2I &h1, TH2I &h2)
3367 hnew.SetDirectory(0);
3375 TH2I operator/(TH2I &h1, TH2I &h2)
3379 hnew.SetDirectory(0);
3395 TH2F::TH2F(): TH2(), TArrayF()
3398 if (fgDefaultSumw2) Sumw2();
3413 TH2F::TH2F(
const char *name,
const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
3414 ,Int_t nbinsy,Double_t ylow,Double_t yup)
3415 :TH2(name,title,nbinsx,xlow,xup,nbinsy,ylow,yup)
3417 TArrayF::Set(fNcells);
3418 if (fgDefaultSumw2) Sumw2();
3420 if (xlow >= xup || ylow >= yup) SetBuffer(fgBufferSize);
3427 TH2F::TH2F(
const char *name,
const char *title,Int_t nbinsx,
const Double_t *xbins
3428 ,Int_t nbinsy,Double_t ylow,Double_t yup)
3429 :TH2(name,title,nbinsx,xbins,nbinsy,ylow,yup)
3431 TArrayF::Set(fNcells);
3432 if (fgDefaultSumw2) Sumw2();
3439 TH2F::TH2F(
const char *name,
const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
3440 ,Int_t nbinsy,
const Double_t *ybins)
3441 :TH2(name,title,nbinsx,xlow,xup,nbinsy,ybins)
3443 TArrayF::Set(fNcells);
3444 if (fgDefaultSumw2) Sumw2();
3451 TH2F::TH2F(
const char *name,
const char *title,Int_t nbinsx,
const Double_t *xbins
3452 ,Int_t nbinsy,
const Double_t *ybins)
3453 :TH2(name,title,nbinsx,xbins,nbinsy,ybins)
3455 TArrayF::Set(fNcells);
3456 if (fgDefaultSumw2) Sumw2();
3463 TH2F::TH2F(
const char *name,
const char *title,Int_t nbinsx,
const Float_t *xbins
3464 ,Int_t nbinsy,
const Float_t *ybins)
3465 :TH2(name,title,nbinsx,xbins,nbinsy,ybins)
3467 TArrayF::Set(fNcells);
3468 if (fgDefaultSumw2) Sumw2();
3475 TH2F::TH2F(
const TMatrixFBase &m)
3476 :TH2(
"TMatrixFBase",
"",m.GetNcols(),m.GetColLwb(),1+m.GetColUpb(),m.GetNrows(),m.GetRowLwb(),1+m.GetRowUpb())
3478 TArrayF::Set(fNcells);
3479 Int_t ilow = m.GetRowLwb();
3480 Int_t iup = m.GetRowUpb();
3481 Int_t jlow = m.GetColLwb();
3482 Int_t jup = m.GetColUpb();
3483 for (Int_t i=ilow;i<=iup;i++) {
3484 for (Int_t j=jlow;j<=jup;j++) {
3485 SetBinContent(j-jlow+1,i-ilow+1,m(i,j));
3494 TH2F::TH2F(
const TH2F &h2f) : TH2(), TArrayF()
3496 ((TH2F&)h2f).Copy(*
this);
3503 void TH2F::Copy(TObject &newth2)
const
3505 TH2::Copy((TH2F&)newth2);
3512 void TH2F::Reset(Option_t *option)
3523 void TH2F::SetBinsLength(Int_t n)
3525 if (n < 0) n = (fXaxis.GetNbins()+2)*(fYaxis.GetNbins()+2);
3534 void TH2F::Streamer(TBuffer &R__b)
3536 if (R__b.IsReading()) {
3538 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
3540 R__b.ReadClassBuffer(TH2F::Class(),
this, R__v, R__s, R__c);
3546 TH1::Streamer(R__b);
3547 TArrayF::Streamer(R__b);
3549 R__b >> fScalefactor;
3554 TH2::Streamer(R__b);
3555 TArrayF::Streamer(R__b);
3556 R__b.CheckByteCount(R__s, R__c, TH2F::IsA());
3561 R__b.WriteClassBuffer(TH2F::Class(),
this);
3569 TH2F& TH2F::operator=(
const TH2F &h1)
3571 if (
this != &h1) ((TH2F&)h1).Copy(*
this);
3579 TH2F operator*(Float_t c1, TH2F &h1)
3583 hnew.SetDirectory(0);
3591 TH2F operator*(TH2F &h1, Float_t c1)
3595 hnew.SetDirectory(0);
3603 TH2F operator+(TH2F &h1, TH2F &h2)
3607 hnew.SetDirectory(0);
3615 TH2F operator-(TH2F &h1, TH2F &h2)
3619 hnew.SetDirectory(0);
3627 TH2F operator*(TH2F &h1, TH2F &h2)
3631 hnew.SetDirectory(0);
3639 TH2F operator/(TH2F &h1, TH2F &h2)
3643 hnew.SetDirectory(0);
3659 TH2D::TH2D(): TH2(), TArrayD()
3662 if (fgDefaultSumw2) Sumw2();
3677 TH2D::TH2D(
const char *name,
const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
3678 ,Int_t nbinsy,Double_t ylow,Double_t yup)
3679 :TH2(name,title,nbinsx,xlow,xup,nbinsy,ylow,yup)
3681 TArrayD::Set(fNcells);
3682 if (fgDefaultSumw2) Sumw2();
3684 if (xlow >= xup || ylow >= yup) SetBuffer(fgBufferSize);
3691 TH2D::TH2D(
const char *name,
const char *title,Int_t nbinsx,
const Double_t *xbins
3692 ,Int_t nbinsy,Double_t ylow,Double_t yup)
3693 :TH2(name,title,nbinsx,xbins,nbinsy,ylow,yup)
3695 TArrayD::Set(fNcells);
3696 if (fgDefaultSumw2) Sumw2();
3703 TH2D::TH2D(
const char *name,
const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
3704 ,Int_t nbinsy,
const Double_t *ybins)
3705 :TH2(name,title,nbinsx,xlow,xup,nbinsy,ybins)
3707 TArrayD::Set(fNcells);
3708 if (fgDefaultSumw2) Sumw2();
3715 TH2D::TH2D(
const char *name,
const char *title,Int_t nbinsx,
const Double_t *xbins
3716 ,Int_t nbinsy,
const Double_t *ybins)
3717 :TH2(name,title,nbinsx,xbins,nbinsy,ybins)
3719 TArrayD::Set(fNcells);
3720 if (fgDefaultSumw2) Sumw2();
3727 TH2D::TH2D(
const char *name,
const char *title,Int_t nbinsx,
const Float_t *xbins
3728 ,Int_t nbinsy,
const Float_t *ybins)
3729 :TH2(name,title,nbinsx,xbins,nbinsy,ybins)
3731 TArrayD::Set(fNcells);
3732 if (fgDefaultSumw2) Sumw2();
3739 TH2D::TH2D(
const TMatrixDBase &m)
3740 :TH2(
"TMatrixDBase",
"",m.GetNcols(),m.GetColLwb(),1+m.GetColUpb(),m.GetNrows(),m.GetRowLwb(),1+m.GetRowUpb())
3742 TArrayD::Set(fNcells);
3743 Int_t ilow = m.GetRowLwb();
3744 Int_t iup = m.GetRowUpb();
3745 Int_t jlow = m.GetColLwb();
3746 Int_t jup = m.GetColUpb();
3747 for (Int_t i=ilow;i<=iup;i++) {
3748 for (Int_t j=jlow;j<=jup;j++) {
3749 SetBinContent(j-jlow+1,i-ilow+1,m(i,j));
3752 if (fgDefaultSumw2) Sumw2();
3759 TH2D::TH2D(
const TH2D &h2d) : TH2(), TArrayD()
3761 ((TH2D&)h2d).Copy(*
this);
3768 void TH2D::Copy(TObject &newth2)
const
3770 TH2::Copy((TH2D&)newth2);
3777 void TH2D::Reset(Option_t *option)
3788 void TH2D::SetBinsLength(Int_t n)
3790 if (n < 0) n = (fXaxis.GetNbins()+2)*(fYaxis.GetNbins()+2);
3799 void TH2D::Streamer(TBuffer &R__b)
3801 if (R__b.IsReading()) {
3803 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
3805 R__b.ReadClassBuffer(TH2D::Class(),
this, R__v, R__s, R__c);
3811 TH1::Streamer(R__b);
3812 TArrayD::Streamer(R__b);
3814 R__b >> fScalefactor;
3819 TH2::Streamer(R__b);
3820 TArrayD::Streamer(R__b);
3821 R__b.CheckByteCount(R__s, R__c, TH2D::IsA());
3826 R__b.WriteClassBuffer(TH2D::Class(),
this);
3834 TH2D& TH2D::operator=(
const TH2D &h1)
3836 if (
this != &h1) ((TH2D&)h1).Copy(*
this);
3845 TH2D operator*(Float_t c1, TH2D &h1)
3849 hnew.SetDirectory(0);
3857 TH2D operator+(TH2D &h1, TH2D &h2)
3861 hnew.SetDirectory(0);
3869 TH2D operator-(TH2D &h1, TH2D &h2)
3873 hnew.SetDirectory(0);
3881 TH2D operator*(TH2D &h1, TH2D &h2)
3885 hnew.SetDirectory(0);
3893 TH2D operator/(TH2D &h1, TH2D &h2)
3897 hnew.SetDirectory(0);