#include "TROOT.h"
#include "TClass.h"
#include "THashList.h"
#include "TH2.h"
#include "TVirtualPad.h"
#include "TF2.h"
#include "TProfile.h"
#include "TRandom.h"
#include "TMatrixFBase.h"
#include "TMatrixDBase.h"
#include "THLimitsFinder.h"
#include "TError.h"
#include "TMath.h"
#include "TObjString.h"
#include "TVirtualHistPainter.h"
ClassImp(TH2)
TH2::TH2()
{
fDimension = 2;
fScalefactor = 1;
fTsumwy = fTsumwy2 = fTsumwxy = 0;
}
TH2::TH2(const char *name,const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
,Int_t nbinsy,Double_t ylow,Double_t yup)
:TH1(name,title,nbinsx,xlow,xup)
{
fDimension = 2;
fScalefactor = 1;
fTsumwy = fTsumwy2 = fTsumwxy = 0;
if (nbinsy <= 0) {Warning("TH2","nbinsy is <=0 - set to nbinsy = 1"); nbinsy = 1; }
fYaxis.Set(nbinsy,ylow,yup);
fNcells = fNcells*(nbinsy+2);
}
TH2::TH2(const char *name,const char *title,Int_t nbinsx,const Double_t *xbins
,Int_t nbinsy,Double_t ylow,Double_t yup)
:TH1(name,title,nbinsx,xbins)
{
fDimension = 2;
fScalefactor = 1;
fTsumwy = fTsumwy2 = fTsumwxy = 0;
if (nbinsy <= 0) {Warning("TH2","nbinsy is <=0 - set to nbinsy = 1"); nbinsy = 1; }
fYaxis.Set(nbinsy,ylow,yup);
fNcells = fNcells*(nbinsy+2);
}
TH2::TH2(const char *name,const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
,Int_t nbinsy,const Double_t *ybins)
:TH1(name,title,nbinsx,xlow,xup)
{
fDimension = 2;
fScalefactor = 1;
fTsumwy = fTsumwy2 = fTsumwxy = 0;
if (nbinsy <= 0) {Warning("TH2","nbinsy is <=0 - set to nbinsy = 1"); nbinsy = 1; }
if (ybins) fYaxis.Set(nbinsy,ybins);
else fYaxis.Set(nbinsy,0,1);
fNcells = fNcells*(nbinsy+2);
}
TH2::TH2(const char *name,const char *title,Int_t nbinsx,const Double_t *xbins
,Int_t nbinsy,const Double_t *ybins)
:TH1(name,title,nbinsx,xbins)
{
fDimension = 2;
fScalefactor = 1;
fTsumwy = fTsumwy2 = fTsumwxy = 0;
if (nbinsy <= 0) {Warning("TH2","nbinsy is <=0 - set to nbinsy = 1"); nbinsy = 1; }
if (ybins) fYaxis.Set(nbinsy,ybins);
else fYaxis.Set(nbinsy,0,1);
fNcells = fNcells*(nbinsy+2);
}
TH2::TH2(const char *name,const char *title,Int_t nbinsx,const Float_t *xbins
,Int_t nbinsy,const Float_t *ybins)
:TH1(name,title,nbinsx,xbins)
{
fDimension = 2;
fScalefactor = 1;
fTsumwy = fTsumwy2 = fTsumwxy = 0;
if (nbinsy <= 0) {Warning("TH2","nbinsy is <=0 - set to nbinsy = 1"); nbinsy = 1; }
if (ybins) fYaxis.Set(nbinsy,ybins);
else fYaxis.Set(nbinsy,0,1);
fNcells = fNcells*(nbinsy+2);
}
TH2::TH2(const TH2 &h) : TH1()
{
((TH2&)h).Copy(*this);
}
TH2::~TH2()
{
}
Int_t TH2::BufferEmpty(Int_t action)
{
if (!fBuffer) return 0;
Int_t nbentries = (Int_t)fBuffer[0];
if (nbentries == 0) return 0;
if (nbentries < 0 && action == 0) return 0;
Double_t *buffer = fBuffer;
if (nbentries < 0) {
nbentries = -nbentries;
fBuffer=0;
Reset("ICES");
fBuffer = buffer;
}
if (TestBit(kCanRebin) || fXaxis.GetXmax() <= fXaxis.GetXmin() || fYaxis.GetXmax() <= fYaxis.GetXmin()) {
Double_t xmin = fBuffer[2];
Double_t xmax = xmin;
Double_t ymin = fBuffer[3];
Double_t ymax = ymin;
for (Int_t i=1;i<nbentries;i++) {
Double_t x = fBuffer[3*i+2];
if (x < xmin) xmin = x;
if (x > xmax) xmax = x;
Double_t y = fBuffer[3*i+3];
if (y < ymin) ymin = y;
if (y > ymax) ymax = y;
}
if (fXaxis.GetXmax() <= fXaxis.GetXmin() || fYaxis.GetXmax() <= fYaxis.GetXmin()) {
THLimitsFinder::GetLimitsFinder()->FindGoodLimits(this,xmin,xmax,ymin,ymax);
} else {
fBuffer = 0;
Int_t keep = fBufferSize; fBufferSize = 0;
if (xmin < fXaxis.GetXmin()) RebinAxis(xmin,&fXaxis);
if (xmax >= fXaxis.GetXmax()) RebinAxis(xmax,&fXaxis);
if (ymin < fYaxis.GetXmin()) RebinAxis(ymin,&fYaxis);
if (ymax >= fYaxis.GetXmax()) RebinAxis(ymax,&fYaxis);
fBuffer = buffer;
fBufferSize = keep;
}
}
fBuffer = 0;
for (Int_t i=0;i<nbentries;i++) {
Fill(buffer[3*i+2],buffer[3*i+3],buffer[3*i+1]);
}
fBuffer = buffer;
if (action > 0) { delete [] fBuffer; fBuffer = 0; fBufferSize = 0;}
else {
if (nbentries == (Int_t)fEntries) fBuffer[0] = -nbentries;
else fBuffer[0] = 0;
}
return nbentries;
}
Int_t TH2::BufferFill(Double_t x, Double_t y, Double_t w)
{
if (!fBuffer) return -3;
Int_t nbentries = (Int_t)fBuffer[0];
if (nbentries < 0) {
nbentries = -nbentries;
fBuffer[0] = nbentries;
if (fEntries > 0) {
Double_t *buffer = fBuffer; fBuffer=0;
Reset("ICES");
fBuffer = buffer;
}
}
if (3*nbentries+3 >= fBufferSize) {
BufferEmpty(1);
return Fill(x,y,w);
}
fBuffer[3*nbentries+1] = w;
fBuffer[3*nbentries+2] = x;
fBuffer[3*nbentries+3] = y;
fBuffer[0] += 1;
return -3;
}
void TH2::Copy(TObject &obj) const
{
TH1::Copy(obj);
((TH2&)obj).fScalefactor = fScalefactor;
((TH2&)obj).fTsumwy = fTsumwy;
((TH2&)obj).fTsumwy2 = fTsumwy2;
((TH2&)obj).fTsumwxy = fTsumwxy;
}
Int_t TH2::Fill(Double_t )
{
Error("Fill", "Invalid signature - do nothing");
return -1;
}
Int_t TH2::Fill(Double_t x,Double_t y)
{
if (fBuffer) return BufferFill(x,y,1);
Int_t binx, biny, bin;
fEntries++;
binx = fXaxis.FindBin(x);
biny = fYaxis.FindBin(y);
if (binx <0 || biny <0) return -1;
bin = biny*(fXaxis.GetNbins()+2) + binx;
AddBinContent(bin);
if (fSumw2.fN) ++fSumw2.fArray[bin];
if (binx == 0 || binx > fXaxis.GetNbins()) {
if (!fgStatOverflows) return -1;
}
if (biny == 0 || biny > fYaxis.GetNbins()) {
if (!fgStatOverflows) return -1;
}
++fTsumw;
++fTsumw2;
fTsumwx += x;
fTsumwx2 += x*x;
fTsumwy += y;
fTsumwy2 += y*y;
fTsumwxy += x*y;
return bin;
}
Int_t TH2::Fill(Double_t x, Double_t y, Double_t w)
{
if (fBuffer) return BufferFill(x,y,w);
Int_t binx, biny, bin;
fEntries++;
binx = fXaxis.FindBin(x);
biny = fYaxis.FindBin(y);
if (binx <0 || biny <0) return -1;
bin = biny*(fXaxis.GetNbins()+2) + binx;
AddBinContent(bin,w);
if (fSumw2.fN) fSumw2.fArray[bin] += w*w;
if (binx == 0 || binx > fXaxis.GetNbins()) {
if (!fgStatOverflows) return -1;
}
if (biny == 0 || biny > fYaxis.GetNbins()) {
if (!fgStatOverflows) return -1;
}
Double_t z= w;
fTsumw += z;
fTsumw2 += z*z;
fTsumwx += z*x;
fTsumwx2 += z*x*x;
fTsumwy += z*y;
fTsumwy2 += z*y*y;
fTsumwxy += z*x*y;
return bin;
}
Int_t TH2::Fill(const char *namex, const char *namey, Double_t w)
{
Int_t binx, biny, bin;
fEntries++;
binx = fXaxis.FindBin(namex);
biny = fYaxis.FindBin(namey);
if (binx <0 || biny <0) return -1;
bin = biny*(fXaxis.GetNbins()+2) + binx;
AddBinContent(bin,w);
if (fSumw2.fN) fSumw2.fArray[bin] += w*w;
if (binx == 0 || binx > fXaxis.GetNbins()) return -1;
if (biny == 0 || biny > fYaxis.GetNbins()) return -1;
Double_t x = fXaxis.GetBinCenter(binx);
Double_t y = fYaxis.GetBinCenter(biny);
Double_t z= w;
fTsumw += z;
fTsumw2 += z*z;
fTsumwx += z*x;
fTsumwx2 += z*x*x;
fTsumwy += z*y;
fTsumwy2 += z*y*y;
fTsumwxy += z*x*y;
return bin;
}
Int_t TH2::Fill(const char *namex, Double_t y, Double_t w)
{
Int_t binx, biny, bin;
fEntries++;
binx = fXaxis.FindBin(namex);
biny = fYaxis.FindBin(y);
if (binx <0 || biny <0) return -1;
bin = biny*(fXaxis.GetNbins()+2) + binx;
AddBinContent(bin,w);
if (fSumw2.fN) fSumw2.fArray[bin] += w*w;
if (binx == 0 || binx > fXaxis.GetNbins()) return -1;
if (biny == 0 || biny > fYaxis.GetNbins()) {
if (!fgStatOverflows) return -1;
}
Double_t x = fXaxis.GetBinCenter(binx);
Double_t z= w;
fTsumw += z;
fTsumw2 += z*z;
fTsumwx += z*x;
fTsumwx2 += z*x*x;
fTsumwy += z*y;
fTsumwy2 += z*y*y;
fTsumwxy += z*x*y;
return bin;
}
Int_t TH2::Fill(Double_t x, const char *namey, Double_t w)
{
Int_t binx, biny, bin;
fEntries++;
binx = fXaxis.FindBin(x);
biny = fYaxis.FindBin(namey);
if (binx <0 || biny <0) return -1;
bin = biny*(fXaxis.GetNbins()+2) + binx;
AddBinContent(bin,w);
if (fSumw2.fN) fSumw2.fArray[bin] += w*w;
if (binx == 0 || binx > fXaxis.GetNbins()) {
if (!fgStatOverflows) return -1;
}
if (biny == 0 || biny > fYaxis.GetNbins()) return -1;
Double_t y = fYaxis.GetBinCenter(biny);
Double_t z= w;
fTsumw += z;
fTsumw2 += z*z;
fTsumwx += z*x;
fTsumwx2 += z*x*x;
fTsumwy += z*y;
fTsumwy2 += z*y*y;
fTsumwxy += z*x*y;
return bin;
}
void TH2::FillN(Int_t ntimes, const Double_t *x, const Double_t *y, const Double_t *w, Int_t stride)
{
Int_t binx, biny, bin, i;
ntimes *= stride;
Int_t ifirst = 0;
if (fBuffer) {
for (i=0;i<ntimes;i+=stride) {
if (!fBuffer) break;
if (w) BufferFill(x[i],y[i],w[i]);
else BufferFill(x[i], y[i], 1.);
}
if (i < ntimes && fBuffer==0)
ifirst = i;
else
return;
}
Double_t ww = 1;
for (i=ifirst;i<ntimes;i+=stride) {
fEntries++;
binx = fXaxis.FindBin(x[i]);
biny = fYaxis.FindBin(y[i]);
if (binx <0 || biny <0) continue;
bin = biny*(fXaxis.GetNbins()+2) + binx;
if (w) ww = w[i];
AddBinContent(bin,ww);
if (fSumw2.fN) fSumw2.fArray[bin] += ww*ww;
if (binx == 0 || binx > fXaxis.GetNbins()) {
if (!fgStatOverflows) continue;
}
if (biny == 0 || biny > fYaxis.GetNbins()) {
if (!fgStatOverflows) continue;
}
Double_t z= ww;
fTsumw += z;
fTsumw2 += z*z;
fTsumwx += z*x[i];
fTsumwx2 += z*x[i]*x[i];
fTsumwy += z*y[i];
fTsumwy2 += z*y[i]*y[i];
fTsumwxy += z*x[i]*y[i];
}
}
void TH2::FillRandom(const char *fname, Int_t ntimes)
{
Int_t bin, binx, biny, ibin, loop;
Double_t r1, x, y;
TObject *fobj = gROOT->GetFunction(fname);
if (!fobj) { Error("FillRandom", "Unknown function: %s",fname); return; }
TF2 * f1 = dynamic_cast<TF2*>(fobj);
if (!f1) { Error("FillRandom", "Function: %s is not a TF2",fname); return; }
Int_t nbinsx = GetNbinsX();
Int_t nbinsy = GetNbinsY();
Int_t nbins = nbinsx*nbinsy;
Double_t *integral = new Double_t[nbins+1];
ibin = 0;
integral[ibin] = 0;
for (biny=1;biny<=nbinsy;biny++) {
for (binx=1;binx<=nbinsx;binx++) {
ibin++;
Double_t fint = f1->Integral(fXaxis.GetBinLowEdge(binx), fXaxis.GetBinUpEdge(binx), fYaxis.GetBinLowEdge(biny), fYaxis.GetBinUpEdge(biny));
integral[ibin] = integral[ibin-1] + fint;
}
}
if (integral[nbins] == 0 ) {
delete [] integral;
Error("FillRandom", "Integral = zero"); return;
}
for (bin=1;bin<=nbins;bin++) integral[bin] /= integral[nbins];
for (loop=0;loop<ntimes;loop++) {
r1 = gRandom->Rndm(loop);
ibin = TMath::BinarySearch(nbins,&integral[0],r1);
biny = ibin/nbinsx;
binx = 1 + ibin - nbinsx*biny;
biny++;
x = fXaxis.GetBinCenter(binx);
y = fYaxis.GetBinCenter(biny);
Fill(x,y, 1.);
}
delete [] integral;
}
void TH2::FillRandom(TH1 *h, Int_t ntimes)
{
if (!h) { Error("FillRandom", "Null histogram"); return; }
if (fDimension != h->GetDimension()) {
Error("FillRandom", "Histograms with different dimensions"); return;
}
if (h->ComputeIntegral() == 0) return;
Int_t loop;
Double_t x,y;
TH2 *h2 = (TH2*)h;
for (loop=0;loop<ntimes;loop++) {
h2->GetRandom2(x,y);
Fill(x,y,1.);
}
}
Int_t TH2::FindFirstBinAbove(Double_t threshold, Int_t axis) const
{
if (axis < 1 || axis > 2) {
Warning("FindFirstBinAbove","Invalid axis number : %d, axis x assumed\n",axis);
axis = 1;
}
Int_t nbinsx = fXaxis.GetNbins();
Int_t nbinsy = fYaxis.GetNbins();
Int_t binx, biny;
if (axis == 1) {
for (binx=1;binx<=nbinsx;binx++) {
for (biny=1;biny<=nbinsy;biny++) {
if (GetBinContent(binx,biny) > threshold) return binx;
}
}
} else {
for (biny=1;biny<=nbinsy;biny++) {
for (binx=1;binx<=nbinsx;binx++) {
if (GetBinContent(binx,biny) > threshold) return biny;
}
}
}
return -1;
}
Int_t TH2::FindLastBinAbove(Double_t threshold, Int_t axis) const
{
if (axis < 1 || axis > 2) {
Warning("FindLastBinAbove","Invalid axis number : %d, axis x assumed\n",axis);
axis = 1;
}
Int_t nbinsx = fXaxis.GetNbins();
Int_t nbinsy = fYaxis.GetNbins();
Int_t binx, biny;
if (axis == 1) {
for (binx=nbinsx;binx>=1;binx--) {
for (biny=1;biny<=nbinsy;biny++) {
if (GetBinContent(binx,biny) > threshold) return binx;
}
}
} else {
for (biny=nbinsy;biny>=1;biny--) {
for (binx=1;binx<=nbinsx;binx++) {
if (GetBinContent(binx,biny) > threshold) return biny;
}
}
}
return -1;
}
void TH2::DoFitSlices(bool onX,
TF1 *f1, Int_t firstbin, Int_t lastbin, Int_t cut, Option_t *option, TObjArray* arr)
{
TAxis& outerAxis = (onX ? fYaxis : fXaxis);
TAxis& innerAxis = (onX ? fXaxis : fYaxis);
Int_t nbins = outerAxis.GetNbins();
if (firstbin < 0) firstbin = 0;
if (lastbin < 0 || lastbin > nbins + 1) lastbin = nbins + 1;
if (lastbin < firstbin) {firstbin = 0; lastbin = nbins + 1;}
TString opt = option;
opt.ToLower();
Int_t ngroup = 1;
if (opt.Contains("g2")) {ngroup = 2; opt.ReplaceAll("g2","");}
if (opt.Contains("g3")) {ngroup = 3; opt.ReplaceAll("g3","");}
if (opt.Contains("g4")) {ngroup = 4; opt.ReplaceAll("g4","");}
if (opt.Contains("g5")) {ngroup = 5; opt.ReplaceAll("g5","");}
Int_t nstep = ngroup;
if (opt.Contains("s")) nstep = 1;
if (f1 == 0) {
f1 = (TF1*)gROOT->GetFunction("gaus");
if (f1 == 0) f1 = new TF1("gaus","gaus",innerAxis.GetXmin(),innerAxis.GetXmax());
else f1->SetRange(innerAxis.GetXmin(),innerAxis.GetXmax());
}
Int_t npar = f1->GetNpar();
if (npar <= 0) return;
Double_t *parsave = new Double_t[npar];
f1->GetParameters(parsave);
if (arr) {
arr->SetOwner();
arr->Expand(npar + 1);
}
Int_t ipar;
TH1D **hlist = new TH1D*[npar];
char *name = new char[2000];
char *title = new char[2000];
const TArrayD *bins = outerAxis.GetXbins();
for (ipar=0;ipar<npar;ipar++) {
snprintf(name,2000,"%s_%d",GetName(),ipar);
snprintf(title,2000,"Fitted value of par[%d]=%s",ipar,f1->GetParName(ipar));
delete gDirectory->FindObject(name);
if (bins->fN == 0) {
hlist[ipar] = new TH1D(name,title, nbins, outerAxis.GetXmin(), outerAxis.GetXmax());
} else {
hlist[ipar] = new TH1D(name,title, nbins,bins->fArray);
}
hlist[ipar]->GetXaxis()->SetTitle(outerAxis.GetTitle());
if (arr)
(*arr)[ipar] = hlist[ipar];
}
snprintf(name,2000,"%s_chi2",GetName());
delete gDirectory->FindObject(name);
TH1D *hchi2 = 0;
if (bins->fN == 0) {
hchi2 = new TH1D(name,"chisquare", nbins, outerAxis.GetXmin(), outerAxis.GetXmax());
} else {
hchi2 = new TH1D(name,"chisquare", nbins, bins->fArray);
}
hchi2->GetXaxis()->SetTitle(outerAxis.GetTitle());
if (arr)
(*arr)[npar] = hchi2;
Int_t bin;
Long64_t nentries;
for (bin=firstbin;bin+ngroup-1<=lastbin;bin += nstep) {
TH1D *hp;
if (onX)
hp= ProjectionX("_temp",bin,bin+ngroup-1,"e");
else
hp= ProjectionY("_temp",bin,bin+ngroup-1,"e");
if (hp == 0) continue;
nentries = Long64_t(hp->GetEntries());
if (nentries == 0 || nentries < cut) {delete hp; continue;}
f1->SetParameters(parsave);
hp->Fit(f1,opt.Data());
Int_t npfits = f1->GetNumberFitPoints();
if (npfits > npar && npfits >= cut) {
Int_t binOn = bin + ngroup/2;
for (ipar=0;ipar<npar;ipar++) {
hlist[ipar]->Fill(outerAxis.GetBinCenter(binOn),f1->GetParameter(ipar));
hlist[ipar]->SetBinError(binOn,f1->GetParError(ipar));
}
hchi2->Fill(outerAxis.GetBinCenter(binOn),f1->GetChisquare()/(npfits-npar));
}
delete hp;
}
delete [] parsave;
delete [] name;
delete [] title;
delete [] hlist;
}
void TH2::FitSlicesX(TF1 *f1, Int_t firstybin, Int_t lastybin, Int_t cut, Option_t *option, TObjArray* arr)
{
DoFitSlices(true, f1, firstybin, lastybin, cut, option, arr);
}
void TH2::FitSlicesY(TF1 *f1, Int_t firstxbin, Int_t lastxbin, Int_t cut, Option_t *option, TObjArray* arr)
{
/*
<img src="gif/fitslicesy.gif">
*/
//End_Html
DoFitSlices(false, f1, firstxbin, lastxbin, cut, option, arr);
}
Double_t TH2::GetBinWithContent2(Double_t c, Int_t &binx, Int_t &biny, Int_t firstxbin, Int_t lastxbin,
Int_t firstybin, Int_t lastybin, Double_t maxdiff) const
{
if (fDimension != 2) {
binx = -1;
biny = -1;
Error("GetBinWithContent2","function is only valid for 2-D histograms");
return 0;
}
if (firstxbin < 0) firstxbin = 1;
if (lastxbin < firstxbin) lastxbin = fXaxis.GetNbins();
if (firstybin < 0) firstybin = 1;
if (lastybin < firstybin) lastybin = fYaxis.GetNbins();
Double_t diff, curmax = 1.e240;
for (Int_t j = firstybin; j <= lastybin; j++) {
for (Int_t i = firstxbin; i <= lastxbin; i++) {
diff = TMath::Abs(GetBinContent(i,j)-c);
if (diff <= 0) {binx = i; biny=j; return diff;}
if (diff < curmax && diff <= maxdiff) {curmax = diff, binx=i; biny=j;}
}
}
return curmax;
}
Double_t TH2::GetCorrelationFactor(Int_t axis1, Int_t axis2) const
{
if (axis1 < 1 || axis2 < 1 || axis1 > 2 || axis2 > 2) {
Error("GetCorrelationFactor","Wrong parameters");
return 0;
}
if (axis1 == axis2) return 1;
Double_t rms1 = GetRMS(axis1);
if (rms1 == 0) return 0;
Double_t rms2 = GetRMS(axis2);
if (rms2 == 0) return 0;
return GetCovariance(axis1,axis2)/rms1/rms2;
}
Double_t TH2::GetCovariance(Int_t axis1, Int_t axis2) const
{
if (axis1 < 1 || axis2 < 1 || axis1 > 2 || axis2 > 2) {
Error("GetCovariance","Wrong parameters");
return 0;
}
Double_t stats[kNstat];
GetStats(stats);
Double_t sumw = stats[0];
Double_t sumwx = stats[2];
Double_t sumwx2 = stats[3];
Double_t sumwy = stats[4];
Double_t sumwy2 = stats[5];
Double_t sumwxy = stats[6];
if (sumw == 0) return 0;
if (axis1 == 1 && axis2 == 1) {
return TMath::Abs(sumwx2/sumw - sumwx/sumw*sumwx/sumw);
}
if (axis1 == 2 && axis2 == 2) {
return TMath::Abs(sumwy2/sumw - sumwy/sumw*sumwy/sumw);
}
return sumwxy/sumw - sumwx/sumw*sumwy/sumw;
}
void TH2::GetRandom2(Double_t &x, Double_t &y)
{
Int_t nbinsx = GetNbinsX();
Int_t nbinsy = GetNbinsY();
Int_t nbins = nbinsx*nbinsy;
Double_t integral;
if (fIntegral) {
if (fIntegral[nbins+1] != fEntries) integral = ComputeIntegral(true);
else integral = fIntegral[nbins];
} else {
integral = ComputeIntegral(true);
}
if (integral == 0 ) { x = 0; y = 0; return;}
if (integral == TMath::QuietNaN() ) { x = TMath::QuietNaN(); y = TMath::QuietNaN(); return;}
Double_t r1 = gRandom->Rndm();
Int_t ibin = TMath::BinarySearch(nbins,fIntegral,(Double_t) r1);
Int_t biny = ibin/nbinsx;
Int_t binx = ibin - nbinsx*biny;
x = fXaxis.GetBinLowEdge(binx+1);
if (r1 > fIntegral[ibin]) x +=
fXaxis.GetBinWidth(binx+1)*(r1-fIntegral[ibin])/(fIntegral[ibin+1] - fIntegral[ibin]);
y = fYaxis.GetBinLowEdge(biny+1) + fYaxis.GetBinWidth(biny+1)*gRandom->Rndm();
}
void TH2::GetStats(Double_t *stats) const
{
if (fBuffer) ((TH2*)this)->BufferEmpty();
Int_t bin, binx, biny;
Double_t w,err;
Double_t x,y;
if ((fTsumw == 0 && fEntries > 0) || fXaxis.TestBit(TAxis::kAxisRange) || fYaxis.TestBit(TAxis::kAxisRange)) {
for (bin=0;bin<7;bin++) stats[bin] = 0;
Int_t firstBinX = fXaxis.GetFirst();
Int_t lastBinX = fXaxis.GetLast();
Int_t firstBinY = fYaxis.GetFirst();
Int_t lastBinY = fYaxis.GetLast();
if (fgStatOverflows) {
if ( !fXaxis.TestBit(TAxis::kAxisRange) ) {
if (firstBinX == 1) firstBinX = 0;
if (lastBinX == fXaxis.GetNbins() ) lastBinX += 1;
}
if ( !fYaxis.TestBit(TAxis::kAxisRange) ) {
if (firstBinY == 1) firstBinY = 0;
if (lastBinY == fYaxis.GetNbins() ) lastBinY += 1;
}
}
for (biny = firstBinY; biny <= lastBinY; biny++) {
y = fYaxis.GetBinCenter(biny);
for (binx = firstBinX; binx <= lastBinX; binx++) {
bin = GetBin(binx,biny);
x = fXaxis.GetBinCenter(binx);
w = GetBinContent(bin);
err = TMath::Abs(GetBinError(bin));
stats[0] += w;
stats[1] += err*err;
stats[2] += w*x;
stats[3] += w*x*x;
stats[4] += w*y;
stats[5] += w*y*y;
stats[6] += w*x*y;
}
}
} else {
stats[0] = fTsumw;
stats[1] = fTsumw2;
stats[2] = fTsumwx;
stats[3] = fTsumwx2;
stats[4] = fTsumwy;
stats[5] = fTsumwy2;
stats[6] = fTsumwxy;
}
}
Double_t TH2::Integral(Option_t *option) const
{
return Integral(fXaxis.GetFirst(),fXaxis.GetLast(),
fYaxis.GetFirst(),fYaxis.GetLast(),option);
}
Double_t TH2::Integral(Int_t firstxbin, Int_t lastxbin, Int_t firstybin, Int_t lastybin, Option_t *option) const
{
double err = 0;
return DoIntegral(firstxbin,lastxbin,firstybin,lastybin,-1,0,err,option);
}
Double_t TH2::IntegralAndError(Int_t firstxbin, Int_t lastxbin, Int_t firstybin, Int_t lastybin, Double_t & error, Option_t *option) const
{
return DoIntegral(firstxbin,lastxbin,firstybin,lastybin,-1,0,error,option,kTRUE);
}
Double_t TH2::Interpolate(Double_t)
{
Error("Interpolate","This function must be called with 2 arguments for a TH2");
return 0;
}
Double_t TH2::Interpolate(Double_t x, Double_t y)
{
Double_t f=0;
Double_t x1=0,x2=0,y1=0,y2=0;
Double_t dx,dy;
Int_t bin_x = fXaxis.FindBin(x);
Int_t bin_y = fYaxis.FindBin(y);
if(bin_x<1 || bin_x>GetNbinsX() || bin_y<1 || bin_y>GetNbinsY()) {
Error("Interpolate","Cannot interpolate outside histogram domain.");
return 0;
}
Int_t quadrant = 0;
dx = fXaxis.GetBinUpEdge(bin_x)-x;
dy = fYaxis.GetBinUpEdge(bin_y)-y;
if (dx<=fXaxis.GetBinWidth(bin_x)/2 && dy<=fYaxis.GetBinWidth(bin_y)/2)
quadrant = 1;
if (dx>fXaxis.GetBinWidth(bin_x)/2 && dy<=fYaxis.GetBinWidth(bin_y)/2)
quadrant = 2;
if (dx>fXaxis.GetBinWidth(bin_x)/2 && dy>fYaxis.GetBinWidth(bin_y)/2)
quadrant = 3;
if (dx<=fXaxis.GetBinWidth(bin_x)/2 && dy>fYaxis.GetBinWidth(bin_y)/2)
quadrant = 4;
switch(quadrant) {
case 1:
x1 = fXaxis.GetBinCenter(bin_x);
y1 = fYaxis.GetBinCenter(bin_y);
x2 = fXaxis.GetBinCenter(bin_x+1);
y2 = fYaxis.GetBinCenter(bin_y+1);
break;
case 2:
x1 = fXaxis.GetBinCenter(bin_x-1);
y1 = fYaxis.GetBinCenter(bin_y);
x2 = fXaxis.GetBinCenter(bin_x);
y2 = fYaxis.GetBinCenter(bin_y+1);
break;
case 3:
x1 = fXaxis.GetBinCenter(bin_x-1);
y1 = fYaxis.GetBinCenter(bin_y-1);
x2 = fXaxis.GetBinCenter(bin_x);
y2 = fYaxis.GetBinCenter(bin_y);
break;
case 4:
x1 = fXaxis.GetBinCenter(bin_x);
y1 = fYaxis.GetBinCenter(bin_y-1);
x2 = fXaxis.GetBinCenter(bin_x+1);
y2 = fYaxis.GetBinCenter(bin_y);
break;
}
Int_t bin_x1 = fXaxis.FindBin(x1);
if(bin_x1<1) bin_x1=1;
Int_t bin_x2 = fXaxis.FindBin(x2);
if(bin_x2>GetNbinsX()) bin_x2=GetNbinsX();
Int_t bin_y1 = fYaxis.FindBin(y1);
if(bin_y1<1) bin_y1=1;
Int_t bin_y2 = fYaxis.FindBin(y2);
if(bin_y2>GetNbinsY()) bin_y2=GetNbinsY();
Int_t bin_q22 = GetBin(bin_x2,bin_y2);
Int_t bin_q12 = GetBin(bin_x1,bin_y2);
Int_t bin_q11 = GetBin(bin_x1,bin_y1);
Int_t bin_q21 = GetBin(bin_x2,bin_y1);
Double_t q11 = GetBinContent(bin_q11);
Double_t q12 = GetBinContent(bin_q12);
Double_t q21 = GetBinContent(bin_q21);
Double_t q22 = GetBinContent(bin_q22);
Double_t d = 1.0*(x2-x1)*(y2-y1);
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);
return f;
}
Double_t TH2::Interpolate(Double_t, Double_t, Double_t)
{
Error("Interpolate","This function must be called with 2 arguments for a TH2");
return 0;
}
Double_t TH2::KolmogorovTest(const TH1 *h2, Option_t *option) const
{
TString opt = option;
opt.ToUpper();
Double_t prb = 0;
TH1 *h1 = (TH1*)this;
if (h2 == 0) return 0;
TAxis *xaxis1 = h1->GetXaxis();
TAxis *xaxis2 = h2->GetXaxis();
TAxis *yaxis1 = h1->GetYaxis();
TAxis *yaxis2 = h2->GetYaxis();
Int_t ncx1 = xaxis1->GetNbins();
Int_t ncx2 = xaxis2->GetNbins();
Int_t ncy1 = yaxis1->GetNbins();
Int_t ncy2 = yaxis2->GetNbins();
if (h1->GetDimension() != 2 || h2->GetDimension() != 2) {
Error("KolmogorovTest","Histograms must be 2-D\n");
return 0;
}
if (ncx1 != ncx2) {
Error("KolmogorovTest","Number of channels in X is different, %d and %d\n",ncx1,ncx2);
return 0;
}
if (ncy1 != ncy2) {
Error("KolmogorovTest","Number of channels in Y is different, %d and %d\n",ncy1,ncy2);
return 0;
}
Bool_t afunc1 = kFALSE;
Bool_t afunc2 = kFALSE;
Double_t difprec = 1e-5;
Double_t diff1 = TMath::Abs(xaxis1->GetXmin() - xaxis2->GetXmin());
Double_t diff2 = TMath::Abs(xaxis1->GetXmax() - xaxis2->GetXmax());
if (diff1 > difprec || diff2 > difprec) {
Error("KolmogorovTest","histograms with different binning along X");
return 0;
}
diff1 = TMath::Abs(yaxis1->GetXmin() - yaxis2->GetXmin());
diff2 = TMath::Abs(yaxis1->GetXmax() - yaxis2->GetXmax());
if (diff1 > difprec || diff2 > difprec) {
Error("KolmogorovTest","histograms with different binning along Y");
return 0;
}
Int_t ibeg = 1, jbeg = 1;
Int_t iend = ncx1, jend = ncy1;
if (opt.Contains("U")) {ibeg = 0; jbeg = 0;}
if (opt.Contains("O")) {iend = ncx1+1; jend = ncy1+1;}
Int_t i,j;
Double_t sum1 = 0;
Double_t sum2 = 0;
Double_t w1 = 0;
Double_t w2 = 0;
for (i = ibeg; i <= iend; i++) {
for (j = jbeg; j <= jend; j++) {
sum1 += h1->GetBinContent(i,j);
sum2 += h2->GetBinContent(i,j);
Double_t ew1 = h1->GetBinError(i,j);
Double_t ew2 = h2->GetBinError(i,j);
w1 += ew1*ew1;
w2 += ew2*ew2;
}
}
if (sum1 == 0) {
Error("KolmogorovTest","Integral is zero for h1=%s\n",h1->GetName());
return 0;
}
if (sum2 == 0) {
Error("KolmogorovTest","Integral is zero for h2=%s\n",h2->GetName());
return 0;
}
Double_t esum1 = 0, esum2 = 0;
if (w1 > 0)
esum1 = sum1 * sum1 / w1;
else
afunc1 = kTRUE;
if (w2 > 0)
esum2 = sum2 * sum2 / w2;
else
afunc2 = kTRUE;
if (afunc2 && afunc1) {
Error("KolmogorovTest","Errors are zero for both histograms\n");
return 0;
}
Double_t s1 = 1/sum1;
Double_t s2 = 1/sum2;
Double_t dfmax1 = 0;
Double_t rsum1=0, rsum2=0;
for (i=ibeg;i<=iend;i++) {
for (j=jbeg;j<=jend;j++) {
rsum1 += s1*h1->GetCellContent(i,j);
rsum2 += s2*h2->GetCellContent(i,j);
dfmax1 = TMath::Max(dfmax1, TMath::Abs(rsum1-rsum2));
}
}
Double_t dfmax2 = 0;
rsum1=0, rsum2=0;
for (j=jbeg;j<=jend;j++) {
for (i=ibeg;i<=iend;i++) {
rsum1 += s1*h1->GetCellContent(i,j);
rsum2 += s2*h2->GetCellContent(i,j);
dfmax2 = TMath::Max(dfmax2, TMath::Abs(rsum1-rsum2));
}
}
Double_t factnm;
if (afunc1) factnm = TMath::Sqrt(esum2);
else if (afunc2) factnm = TMath::Sqrt(esum1);
else factnm = TMath::Sqrt(esum1*sum2/(esum1+esum2));
Double_t dfmax = 0.5*(dfmax1+dfmax2);
Double_t z = dfmax*factnm;
prb = TMath::KolmogorovProb(z);
Double_t prb1 = 0, prb2 = 0;
if (opt.Contains("N") && !(afunc1 || afunc2 ) ) {
prb1 = prb;
Double_t d12 = esum1-esum2;
Double_t chi2 = d12*d12/(esum1+esum2);
prb2 = TMath::Prob(chi2,1);
if (prb > 0 && prb2 > 0) prb = prb*prb2*(1-TMath::Log(prb*prb2));
else prb = 0;
}
if (opt.Contains("D")) {
printf(" Kolmo Prob h1 = %s, sum1=%g\n",h1->GetName(),sum1);
printf(" Kolmo Prob h2 = %s, sum2=%g\n",h2->GetName(),sum2);
printf(" Kolmo Probabil = %f, Max Dist = %g\n",prb,dfmax);
if (opt.Contains("N"))
printf(" Kolmo Probabil = %f for shape alone, =%f for normalisation alone\n",prb1,prb2);
}
if (TMath::Abs(rsum1-1) > 0.002) Warning("KolmogorovTest","Numerical problems with h1=%s\n",h1->GetName());
if (TMath::Abs(rsum2-1) > 0.002) Warning("KolmogorovTest","Numerical problems with h2=%s\n",h2->GetName());
if(opt.Contains("M")) return dfmax;
return prb;
}
Long64_t TH2::Merge(TCollection *list)
{
if (!list) return 0;
if (list->IsEmpty()) return (Long64_t) GetEntries();
TList inlist;
inlist.AddAll(list);
TAxis newXAxis;
TAxis newYAxis;
Bool_t initialLimitsFound = kFALSE;
Bool_t allSameLimits = kTRUE;
Bool_t sameLimitsX = kTRUE;
Bool_t sameLimitsY = kTRUE;
Bool_t allHaveLimits = kTRUE;
Bool_t firstHistWithLimits = kTRUE;
TIter next(&inlist);
TH2 * h = this;
do {
Bool_t hasLimits = h->GetXaxis()->GetXmin() < h->GetXaxis()->GetXmax();
allHaveLimits = allHaveLimits && hasLimits;
if (hasLimits) {
h->BufferEmpty();
if (firstHistWithLimits ) {
if (h != this ) {
if (!SameLimitsAndNBins(fXaxis, *(h->GetXaxis())) ) {
if (h->GetXaxis()->GetXbins()->GetSize() != 0) fXaxis.Set(h->GetXaxis()->GetNbins(), h->GetXaxis()->GetXbins()->GetArray());
else fXaxis.Set(h->GetXaxis()->GetNbins(), h->GetXaxis()->GetXmin(), h->GetXaxis()->GetXmax());
}
if (!SameLimitsAndNBins(fYaxis, *(h->GetYaxis())) ) {
if (h->GetYaxis()->GetXbins()->GetSize() != 0) fYaxis.Set(h->GetYaxis()->GetNbins(), h->GetYaxis()->GetXbins()->GetArray());
else fYaxis.Set(h->GetYaxis()->GetNbins(), h->GetYaxis()->GetXmin(), h->GetYaxis()->GetXmax());
}
}
firstHistWithLimits = kFALSE;
}
if (!initialLimitsFound) {
initialLimitsFound = kTRUE;
if (h->GetXaxis()->GetXbins()->GetSize() != 0) newXAxis.Set(h->GetXaxis()->GetNbins(), h->GetXaxis()->GetXbins()->GetArray());
else newXAxis.Set(h->GetXaxis()->GetNbins(), h->GetXaxis()->GetXmin(), h->GetXaxis()->GetXmax());
if (h->GetYaxis()->GetXbins()->GetSize() != 0) newYAxis.Set(h->GetYaxis()->GetNbins(), h->GetYaxis()->GetXbins()->GetArray());
else newYAxis.Set(h->GetYaxis()->GetNbins(), h->GetYaxis()->GetXmin(), h->GetYaxis()->GetXmax());
}
else {
if (!SameLimitsAndNBins(newXAxis, *(h->GetXaxis()))) {
sameLimitsX = kFALSE;
if (!RecomputeAxisLimits(newXAxis, *(h->GetXaxis()))) {
Error("Merge", "Cannot merge histograms - limits are inconsistent:\n "
"first: (%d, %f, %f), second: (%d, %f, %f)",
newXAxis.GetNbins(), newXAxis.GetXmin(), newXAxis.GetXmax(),
h->GetXaxis()->GetNbins(), h->GetXaxis()->GetXmin(),
h->GetXaxis()->GetXmax());
return -1;
}
}
if (!SameLimitsAndNBins(newYAxis, *(h->GetYaxis()))) {
sameLimitsY = kFALSE;
if (!RecomputeAxisLimits(newYAxis, *(h->GetYaxis()))) {
Error("Merge", "Cannot merge histograms - limits are inconsistent:\n "
"first: (%d, %f, %f), second: (%d, %f, %f)",
newYAxis.GetNbins(), newYAxis.GetXmin(), newYAxis.GetXmax(),
h->GetYaxis()->GetNbins(), h->GetYaxis()->GetXmin(),
h->GetYaxis()->GetXmax());
return -1;
}
}
allSameLimits = sameLimitsY && sameLimitsX;
}
}
} while ( ( h = dynamic_cast<TH2*> ( next() ) ) != NULL );
if (!h && (*next) ) {
Error("Merge","Attempt to merge object of class: %s to a %s",
(*next)->ClassName(),this->ClassName());
return -1;
}
next.Reset();
TH2 * hclone = 0;
if (!allSameLimits) {
Bool_t mustCleanup = TestBit(kMustCleanup);
if (mustCleanup) ResetBit(kMustCleanup);
hclone = (TH2*)IsA()->New();
hclone->SetDirectory(0);
Copy(*hclone);
if (mustCleanup) SetBit(kMustCleanup);
BufferEmpty(1);
Reset();
SetEntries(0);
inlist.AddFirst(hclone);
}
if (!allSameLimits && initialLimitsFound) {
if (!sameLimitsX) {
fXaxis.SetRange(0,0);
if (newXAxis.GetXbins()->GetSize() != 0) fXaxis.Set(newXAxis.GetNbins(),newXAxis.GetXbins()->GetArray());
else fXaxis.Set(newXAxis.GetNbins(),newXAxis.GetXmin(), newXAxis.GetXmax());
}
if (!sameLimitsY) {
fYaxis.SetRange(0,0);
if (newYAxis.GetXbins()->GetSize() != 0) fYaxis.Set(newYAxis.GetNbins(),newYAxis.GetXbins()->GetArray());
else fYaxis.Set(newYAxis.GetNbins(),newYAxis.GetXmin(), newYAxis.GetXmax());
}
fZaxis.Set(1,0,1);
fNcells = (fXaxis.GetNbins()+2)*(fYaxis.GetNbins()+2);
SetBinsLength(fNcells);
if (fSumw2.fN) {
fSumw2.Set(fNcells);
}
}
if (!allHaveLimits) {
while ( (h = dynamic_cast<TH2*> (next())) ) {
if (h->GetXaxis()->GetXmin() >= h->GetXaxis()->GetXmax() && h->fBuffer) {
Int_t nbentries = (Int_t)h->fBuffer[0];
for (Int_t i = 0; i < nbentries; i++)
Fill(h->fBuffer[3*i + 2], h->fBuffer[3*i + 3], h->fBuffer[3*i + 1]);
}
}
if (!initialLimitsFound) {
if (hclone) {
inlist.Remove(hclone);
delete hclone;
}
return (Long64_t) GetEntries();
}
next.Reset();
}
Double_t stats[kNstat], totstats[kNstat];
for (Int_t i=0;i<kNstat;i++) {totstats[i] = stats[i] = 0;}
GetStats(totstats);
Double_t nentries = GetEntries();
Int_t binx, biny, ix, iy, nx, ny, bin, ibin;
Double_t cu;
Int_t nbix = fXaxis.GetNbins();
Bool_t canRebin=TestBit(kCanRebin);
ResetBit(kCanRebin);
while ((h=(TH2*)next())) {
Double_t histEntries = h->GetEntries();
if (h->fTsumw == 0 && histEntries == 0) continue;
if (h->GetXaxis()->GetXmin() < h->GetXaxis()->GetXmax()) {
h->GetStats(stats);
for (Int_t i = 0; i < kNstat; i++)
totstats[i] += stats[i];
nentries += histEntries;
nx = h->GetXaxis()->GetNbins();
ny = h->GetYaxis()->GetNbins();
for (biny = 0; biny <= ny + 1; biny++) {
if (!allSameLimits)
iy = fYaxis.FindBin(h->GetYaxis()->GetBinCenter(biny));
else
iy = biny;
for (binx = 0; binx <= nx + 1; binx++) {
bin = binx +(nx+2)*biny;
cu = h->GetBinContent(bin);
if (!allSameLimits) {
if (cu != 0 && ( (!sameLimitsX && (binx == 0 || binx == nx+1)) || (!sameLimitsY && (biny == 0 || biny == ny+1)) )) {
Error("Merge", "Cannot merge histograms - the histograms have"
" different limits and undeflows/overflows are present."
" The initial histogram is now broken!");
return -1;
}
ix = fXaxis.FindBin(h->GetXaxis()->GetBinCenter(binx));
}
else {
ix = binx;
}
ibin = ix +(nbix+2)*iy;
if (ibin < 0) continue;
AddBinContent(ibin,cu);
if (fSumw2.fN) {
Double_t error1 = h->GetBinError(bin);
fSumw2.fArray[ibin] += error1*error1;
}
}
}
}
}
if (canRebin) SetBit(kCanRebin);
PutStats(totstats);
SetEntries(nentries);
if (hclone) {
inlist.Remove(hclone);
delete hclone;
}
return (Long64_t)nentries;
}
TH2 *TH2::RebinX(Int_t ngroup, const char *newname)
{
return Rebin2D(ngroup, 1, newname);
}
TH2 *TH2::RebinY(Int_t ngroup, const char *newname)
{
return Rebin2D(1, ngroup, newname);
}
TH2 *TH2::Rebin2D(Int_t nxgroup, Int_t nygroup, const char *newname)
{
Int_t i,j,xbin,ybin;
Int_t nxbins = fXaxis.GetNbins();
Int_t nybins = fYaxis.GetNbins();
Double_t xmin = fXaxis.GetXmin();
Double_t xmax = fXaxis.GetXmax();
Double_t ymin = fYaxis.GetXmin();
Double_t ymax = fYaxis.GetXmax();
if ((nxgroup <= 0) || (nxgroup > nxbins)) {
Error("Rebin", "Illegal value of nxgroup=%d",nxgroup);
return 0;
}
if ((nygroup <= 0) || (nygroup > nybins)) {
Error("Rebin", "Illegal value of nygroup=%d",nygroup);
return 0;
}
Int_t newxbins = nxbins/nxgroup;
Int_t newybins = nybins/nygroup;
Double_t entries = fEntries;
Double_t *oldBins = new Double_t[(nxbins+2)*(nybins+2)];
for (xbin = 0; xbin < nxbins+2; xbin++) {
for (ybin = 0; ybin < nybins+2; ybin++) {
oldBins[ybin*(nxbins+2)+xbin] = GetBinContent(xbin, ybin);
}
}
Double_t *oldErrors = 0;
if (fSumw2.fN != 0) {
oldErrors = new Double_t[(nxbins+2)*(nybins+2)];
for (xbin = 0; xbin < nxbins+2; xbin++) {
for (ybin = 0; ybin < nybins+2; ybin++) {
oldErrors[ybin*(nxbins+2)+xbin] = GetBinError(xbin, ybin);
}
}
}
TH2 *hnew = this;
if (newname && strlen(newname)) {
hnew = (TH2*)Clone();
hnew->SetName(newname);
}
Int_t bitRebin = hnew->TestBit(kCanRebin);
hnew->SetBit(kCanRebin,0);
Double_t stat[kNstat];
GetStats(stat);
bool resetStat = false;
if(newxbins*nxgroup != nxbins) {
xmax = fXaxis.GetBinUpEdge(newxbins*nxgroup);
resetStat = true;
}
if(newybins*nygroup != nybins) {
ymax = fYaxis.GetBinUpEdge(newybins*nygroup);
resetStat = true;
}
Int_t nXdivisions = fXaxis.GetNdivisions();
Color_t xAxisColor = fXaxis.GetAxisColor();
Color_t xLabelColor = fXaxis.GetLabelColor();
Style_t xLabelFont = fXaxis.GetLabelFont();
Float_t xLabelOffset = fXaxis.GetLabelOffset();
Float_t xLabelSize = fXaxis.GetLabelSize();
Float_t xTickLength = fXaxis.GetTickLength();
Float_t xTitleOffset = fXaxis.GetTitleOffset();
Float_t xTitleSize = fXaxis.GetTitleSize();
Color_t xTitleColor = fXaxis.GetTitleColor();
Style_t xTitleFont = fXaxis.GetTitleFont();
Int_t nYdivisions = fYaxis.GetNdivisions();
Color_t yAxisColor = fYaxis.GetAxisColor();
Color_t yLabelColor = fYaxis.GetLabelColor();
Style_t yLabelFont = fYaxis.GetLabelFont();
Float_t yLabelOffset = fYaxis.GetLabelOffset();
Float_t yLabelSize = fYaxis.GetLabelSize();
Float_t yTickLength = fYaxis.GetTickLength();
Float_t yTitleOffset = fYaxis.GetTitleOffset();
Float_t yTitleSize = fYaxis.GetTitleSize();
Color_t yTitleColor = fYaxis.GetTitleColor();
Style_t yTitleFont = fYaxis.GetTitleFont();
if (nxgroup != 1 || nygroup != 1) {
if(fXaxis.GetXbins()->GetSize() > 0 || fYaxis.GetXbins()->GetSize() > 0){
Double_t *xbins = new Double_t[newxbins+1];
for(i = 0; i <= newxbins; ++i) xbins[i] = fXaxis.GetBinLowEdge(1+i*nxgroup);
Double_t *ybins = new Double_t[newybins+1];
for(i = 0; i <= newybins; ++i) ybins[i] = fYaxis.GetBinLowEdge(1+i*nygroup);
hnew->SetBins(newxbins,xbins, newybins, ybins);
delete [] xbins;
delete [] ybins;
} else {
hnew->SetBins(newxbins, xmin, xmax, newybins, ymin, ymax);
}
Double_t binContent, binError;
Int_t oldxbin = 1;
Int_t oldybin = 1;
Int_t bin;
for (xbin = 1; xbin <= newxbins; xbin++) {
oldybin = 1;
for (ybin = 1; ybin <= newybins; ybin++) {
binContent = 0;
binError = 0;
for (i = 0; i < nxgroup; i++) {
if (oldxbin+i > nxbins) break;
for (j =0; j < nygroup; j++) {
if (oldybin+j > nybins) break;
bin = oldxbin + i + (oldybin + j)*(nxbins + 2);
binContent += oldBins[bin];
if (oldErrors) binError += oldErrors[bin]*oldErrors[bin];
}
}
hnew->SetBinContent(xbin,ybin, binContent);
if (oldErrors) hnew->SetBinError(xbin,ybin,TMath::Sqrt(binError));
oldybin += nygroup;
}
oldxbin += nxgroup;
}
hnew->SetBinContent(0,0,oldBins[0]);
if (oldErrors) hnew->SetBinError(0,0,oldErrors[0]);
binContent = 0;
binError = 0;
for(xbin = oldxbin; xbin <= nxbins+1; xbin++)
for(ybin = oldybin; ybin <= nybins+1; ybin++){
bin = xbin + (nxbins+2)*ybin;
binContent += oldBins[bin];
if(oldErrors) binError += oldErrors[bin]*oldErrors[bin];
}
hnew->SetBinContent(newxbins+1,newybins+1,binContent);
if(oldErrors) hnew->SetBinError(newxbins+1,newybins+1,TMath::Sqrt(binError));
binContent = 0;
binError = 0;
for(ybin = oldybin; ybin <= nybins+1; ybin++){
bin = ybin*(nxbins+2);
binContent += oldBins[bin];
if(oldErrors) binError += oldErrors[bin] * oldErrors[bin];
}
hnew->SetBinContent(0,newybins+1,binContent);
if(oldErrors) hnew->SetBinError(0,newybins+1,TMath::Sqrt(binError));
binContent = 0;
binError = 0;
for(xbin = oldxbin; xbin <= nxbins+1; xbin++){
bin = xbin;
binContent += oldBins[bin];
if(oldErrors) binError += oldErrors[bin] * oldErrors[bin];
}
hnew->SetBinContent(newxbins+1,0,binContent);
if(oldErrors) hnew->SetBinError(newxbins+1,0,TMath::Sqrt(binError));
Double_t binContent0, binContent2;
Double_t binError0, binError2;
Int_t oldxbin2, oldybin2;
Int_t ufbin, ofbin;
oldxbin2 = 1;
for (xbin = 1; xbin<=newxbins; xbin++) {
binContent0 = binContent2 = 0;
binError0 = binError2 = 0;
for (i=0; i<nxgroup; i++) {
if (oldxbin2+i > nxbins) break;
ufbin = oldxbin2 + i;
binContent0 += oldBins[ufbin];
if(oldErrors) binError0 += oldErrors[ufbin] * oldErrors[ufbin];
for(ybin = oldybin; ybin <= nybins + 1; ybin++){
ofbin = ufbin + ybin*(nxbins+2);
binContent2 += oldBins[ofbin];
if(oldErrors) binError2 += oldErrors[ofbin] * oldErrors[ofbin];
}
}
hnew->SetBinContent(xbin,0,binContent0);
hnew->SetBinContent(xbin,newybins+1,binContent2);
if (oldErrors) {
hnew->SetBinError(xbin,0,TMath::Sqrt(binError0));
hnew->SetBinError(xbin,newybins+1,TMath::Sqrt(binError2) );
}
oldxbin2 += nxgroup;
}
oldybin2 = 1;
for (ybin = 1; ybin<=newybins; ybin++){
binContent0 = binContent2 = 0;
binError0 = binError2 = 0;
for (i=0; i<nygroup; i++) {
if (oldybin2+i > nybins) break;
ufbin = (oldybin2 + i)*(nxbins+2);
binContent0 += oldBins[ufbin];
if(oldErrors) binError0 += oldErrors[ufbin] * oldErrors[ufbin];
for(xbin = oldxbin; xbin <= nxbins+1; xbin++){
ofbin = ufbin + xbin;
binContent2 += oldBins[ofbin];
if(oldErrors) binError2 += oldErrors[ofbin] * oldErrors[ofbin];
}
}
hnew->SetBinContent(0,ybin,binContent0);
hnew->SetBinContent(newxbins+1,ybin,binContent2);
if (oldErrors) {
hnew->SetBinError(0,ybin, TMath::Sqrt(binError0));
hnew->SetBinError(newxbins+1, ybin, TMath::Sqrt(binError2));
}
oldybin2 += nygroup;
}
}
fXaxis.SetNdivisions(nXdivisions);
fXaxis.SetAxisColor(xAxisColor);
fXaxis.SetLabelColor(xLabelColor);
fXaxis.SetLabelFont(xLabelFont);
fXaxis.SetLabelOffset(xLabelOffset);
fXaxis.SetLabelSize(xLabelSize);
fXaxis.SetTickLength(xTickLength);
fXaxis.SetTitleOffset(xTitleOffset);
fXaxis.SetTitleSize(xTitleSize);
fXaxis.SetTitleColor(xTitleColor);
fXaxis.SetTitleFont(xTitleFont);
fYaxis.SetNdivisions(nYdivisions);
fYaxis.SetAxisColor(yAxisColor);
fYaxis.SetLabelColor(yLabelColor);
fYaxis.SetLabelFont(yLabelFont);
fYaxis.SetLabelOffset(yLabelOffset);
fYaxis.SetLabelSize(yLabelSize);
fYaxis.SetTickLength(yTickLength);
fYaxis.SetTitleOffset(yTitleOffset);
fYaxis.SetTitleSize(yTitleSize);
fYaxis.SetTitleColor(yTitleColor);
fYaxis.SetTitleFont(yTitleFont);
hnew->SetEntries(entries);
if (!resetStat) hnew->PutStats(stat);
hnew->SetBit(kCanRebin,bitRebin);
delete [] oldBins;
if (oldErrors) delete [] oldErrors;
return hnew;
}
TProfile *TH2::DoProfile(bool onX, const char *name, Int_t firstbin, Int_t lastbin, Option_t *option) const
{
TString opt = option;
TString cut;
Int_t i1 = opt.Index("[");
if (i1>=0) {
Int_t i2 = opt.Index("]");
cut = opt(i1,i2-i1+1);
}
opt.ToLower();
bool originalRange = opt.Contains("o");
const TAxis& outAxis = ( onX ? fXaxis : fYaxis );
const TAxis& inAxis = ( onX ? fYaxis : fXaxis );
Int_t inN = inAxis.GetNbins();
const char *expectedName = ( onX ? "_pfx" : "_pfy" );
Int_t firstOutBin, lastOutBin;
firstOutBin = outAxis.GetFirst();
lastOutBin = outAxis.GetLast();
if (firstOutBin == 0 && lastOutBin == 0) {
firstOutBin = 1; lastOutBin = outAxis.GetNbins();
}
if ( lastbin < firstbin && inAxis.TestBit(TAxis::kAxisRange) ) {
firstbin = inAxis.GetFirst();
lastbin = inAxis.GetLast();
if (firstbin == 0 && lastbin == 0)
{
firstbin = 1;
lastbin = inAxis.GetNbins();
}
}
if (firstbin < 0) firstbin = 1;
if (lastbin < 0) lastbin = inN;
if (lastbin > inN+1) lastbin = inN;
char *pname = (char*)name;
if (name && strcmp(name, expectedName) == 0) {
Int_t nch = strlen(GetName()) + 5;
pname = new char[nch];
snprintf(pname,nch,"%s%s",GetName(),name);
}
TProfile *h1=0;
TObject *h1obj = gROOT->FindObject(pname);
if (h1obj && h1obj->InheritsFrom(TH1::Class())) {
if (h1obj->IsA() != TProfile::Class() ) {
Error("DoProfile","Histogram with name %s must be a TProfile and is a %s",name,h1obj->ClassName());
return 0;
}
h1 = (TProfile*)h1obj;
h1->Reset();
const TArrayD *xbins = outAxis.GetXbins();
if (xbins->fN == 0) {
if ( originalRange )
h1->SetBins(outAxis.GetNbins(),outAxis.GetXmin(),outAxis.GetXmax());
else
h1->SetBins(lastOutBin-firstOutBin+1,outAxis.GetBinLowEdge(firstOutBin),outAxis.GetBinUpEdge(lastOutBin));
} else {
if (originalRange )
h1->SetBins(outAxis.GetNbins(),xbins->fArray);
else
h1->SetBins(lastOutBin-firstOutBin+1,&xbins->fArray[firstOutBin-1]);
}
}
Int_t ncuts = 0;
if (opt.Contains("[")) {
((TH2 *)this)->GetPainter();
if (fPainter) ncuts = fPainter->MakeCuts((char*)cut.Data());
}
if (!h1) {
const TArrayD *bins = outAxis.GetXbins();
if (bins->fN == 0) {
if ( originalRange )
h1 = new TProfile(pname,GetTitle(),outAxis.GetNbins(),outAxis.GetXmin(),outAxis.GetXmax(),opt);
else
h1 = new TProfile(pname,GetTitle(),lastOutBin-firstOutBin+1,
outAxis.GetBinLowEdge(firstOutBin),
outAxis.GetBinUpEdge(lastOutBin), opt);
} else {
if (originalRange )
h1 = new TProfile(pname,GetTitle(),outAxis.GetNbins(),bins->fArray,opt);
else
h1 = new TProfile(pname,GetTitle(),lastOutBin-firstOutBin+1,&bins->fArray[firstOutBin-1],opt);
}
}
if (pname != name) delete [] pname;
h1->GetXaxis()->ImportAttributes( &outAxis);
h1->SetLineColor(this->GetLineColor());
h1->SetFillColor(this->GetFillColor());
h1->SetMarkerColor(this->GetMarkerColor());
h1->SetMarkerStyle(this->GetMarkerStyle());
bool useWeights = (GetSumw2N() > 0);
if (useWeights) h1->Sumw2();
Double_t totcont = 0;
TArrayD & binSumw2 = *(h1->GetBinSumw2());
for ( Int_t outbin = 0; outbin <= outAxis.GetNbins() + 1; ++outbin) {
if (outAxis.TestBit(TAxis::kAxisRange) && ( outbin < firstOutBin || outbin > lastOutBin )) continue;
Double_t xOut = outAxis.GetBinCenter(outbin);
Int_t binOut = h1->GetXaxis()->FindBin( xOut );
if (binOut <0) continue;
for (Int_t inbin = firstbin ; inbin <= lastbin ; ++inbin) {
Int_t binx, biny;
if (onX) { binx = outbin; biny=inbin; }
else { binx = inbin; biny=outbin; }
if (ncuts) {
if (!fPainter->IsInside(binx,biny)) continue;
}
Int_t bin = GetBin(binx, biny);
Double_t cxy = GetBinContent(bin);
if (cxy) {
Double_t tmp = 0;
if ( useWeights ) tmp = binSumw2.fArray[binOut];
h1->Fill( xOut, inAxis.GetBinCenter(inbin), cxy );
if ( useWeights ) binSumw2.fArray[binOut] = tmp + fSumw2.fArray[bin];
totcont += cxy;
}
}
}
h1->ResetStats();
h1->SetEntries( h1->GetEffectiveEntries() );
if (opt.Contains("d")) {
TVirtualPad *padsav = gPad;
TVirtualPad *pad = gROOT->GetSelectedPad();
if (pad) pad->cd();
opt.Remove(opt.First("d"),1);
if (!gPad || !gPad->FindObject(h1)) {
h1->Draw(opt);
} else {
h1->Paint(opt);
}
if (padsav) padsav->cd();
}
return h1;
}
TProfile *TH2::ProfileX(const char *name, Int_t firstybin, Int_t lastybin, Option_t *option) const
{
return DoProfile(true, name, firstybin, lastybin, option);
}
TProfile *TH2::ProfileY(const char *name, Int_t firstxbin, Int_t lastxbin, Option_t *option) const
{
return DoProfile(false, name, firstxbin, lastxbin, option);
}
TH1D *TH2::DoProjection(bool onX, const char *name, Int_t firstbin, Int_t lastbin, Option_t *option) const
{
const char *expectedName = 0;
Int_t inNbin;
Int_t firstOutBin, lastOutBin;
TAxis* outAxis;
TAxis* inAxis;
TString opt = option;
TString cut;
Int_t i1 = opt.Index("[");
if (i1>=0) {
Int_t i2 = opt.Index("]");
cut = opt(i1,i2-i1+1);
}
opt.ToLower();
bool originalRange = opt.Contains("o");
if ( onX )
{
expectedName = "_px";
inNbin = fYaxis.GetNbins();
outAxis = GetXaxis();
inAxis = GetYaxis();
}
else
{
expectedName = "_py";
inNbin = fXaxis.GetNbins();
outAxis = GetYaxis();
inAxis = GetXaxis();
}
firstOutBin = outAxis->GetFirst();
lastOutBin = outAxis->GetLast();
if (firstOutBin == 0 && lastOutBin == 0) {
firstOutBin = 1; lastOutBin = outAxis->GetNbins();
}
if ( lastbin < firstbin && inAxis->TestBit(TAxis::kAxisRange) ) {
firstbin = inAxis->GetFirst();
lastbin = inAxis->GetLast();
if (firstbin == 0 && lastbin == 0)
{
firstbin = 1;
lastbin = inAxis->GetNbins();
}
}
if (firstbin < 0) firstbin = 0;
if (lastbin < 0) lastbin = inNbin + 1;
if (lastbin > inNbin+1) lastbin = inNbin + 1;
char *pname = (char*)name;
if (name && strcmp(name,expectedName) == 0) {
Int_t nch = strlen(GetName()) + 4;
pname = new char[nch];
snprintf(pname,nch,"%s%s",GetName(),name);
}
TH1D *h1=0;
TObject *h1obj = gROOT->FindObject(pname);
if (h1obj && h1obj->InheritsFrom(TH1::Class())) {
if (h1obj->IsA() != TH1D::Class() ) {
Error("DoProjection","Histogram with name %s must be a TH1D and is a %s",name,h1obj->ClassName());
return 0;
}
h1 = (TH1D*)h1obj;
h1->Reset();
const TArrayD *xbins = outAxis->GetXbins();
if (xbins->fN == 0) {
if ( originalRange )
h1->SetBins(outAxis->GetNbins(),outAxis->GetXmin(),outAxis->GetXmax());
else
h1->SetBins(lastOutBin-firstOutBin+1,outAxis->GetBinLowEdge(firstOutBin),outAxis->GetBinUpEdge(lastOutBin));
} else {
if (originalRange )
h1->SetBins(outAxis->GetNbins(),xbins->fArray);
else
h1->SetBins(lastOutBin-firstOutBin+1,&xbins->fArray[firstOutBin-1]);
}
}
Int_t ncuts = 0;
if (opt.Contains("[")) {
((TH2 *)this)->GetPainter();
if (fPainter) ncuts = fPainter->MakeCuts((char*)cut.Data());
}
if (!h1) {
const TArrayD *bins = outAxis->GetXbins();
if (bins->fN == 0) {
if ( originalRange )
h1 = new TH1D(pname,GetTitle(),outAxis->GetNbins(),outAxis->GetXmin(),outAxis->GetXmax());
else
h1 = new TH1D(pname,GetTitle(),lastOutBin-firstOutBin+1,
outAxis->GetBinLowEdge(firstOutBin),outAxis->GetBinUpEdge(lastOutBin));
} else {
if (originalRange )
h1 = new TH1D(pname,GetTitle(),outAxis->GetNbins(),bins->fArray);
else
h1 = new TH1D(pname,GetTitle(),lastOutBin-firstOutBin+1,&bins->fArray[firstOutBin-1]);
}
if (opt.Contains("e") || GetSumw2N() ) h1->Sumw2();
}
if (pname != name) delete [] pname;
h1->GetXaxis()->ImportAttributes(outAxis);
THashList* labels=outAxis->GetLabels();
if (labels) {
TIter iL(labels);
TObjString* lb;
Int_t i = 1;
while ((lb=(TObjString*)iL())) {
h1->GetXaxis()->SetBinLabel(i,lb->String().Data());
i++;
}
}
h1->SetLineColor(this->GetLineColor());
h1->SetFillColor(this->GetFillColor());
h1->SetMarkerColor(this->GetMarkerColor());
h1->SetMarkerStyle(this->GetMarkerStyle());
Double_t cont,err2;
Double_t totcont = 0;
Bool_t computeErrors = h1->GetSumw2N();
for ( Int_t outbin = 0; outbin <= outAxis->GetNbins() + 1; ++outbin) {
err2 = 0;
cont = 0;
if (outAxis->TestBit(TAxis::kAxisRange) && ( outbin < firstOutBin || outbin > lastOutBin )) continue;
for (Int_t inbin = firstbin ; inbin <= lastbin ; ++inbin) {
Int_t binx, biny;
if (onX) { binx = outbin; biny=inbin; }
else { binx = inbin; biny=outbin; }
if (ncuts) {
if (!fPainter->IsInside(binx,biny)) continue;
}
cont += GetCellContent(binx,biny);
if (computeErrors) {
Double_t exy = GetCellError(binx,biny);
err2 += exy*exy;
}
}
Int_t binOut = h1->GetXaxis()->FindBin( outAxis->GetBinCenter(outbin) );
h1->SetBinContent(binOut ,cont);
if (computeErrors) h1->SetBinError(binOut,TMath::Sqrt(err2));
totcont += cont;
}
bool reuseStats = false;
if ( ( fgStatOverflows == false && firstbin == 1 && lastbin == inNbin ) ||
( fgStatOverflows == true && firstbin == 0 && lastbin == inNbin + 1 ) )
reuseStats = true;
else {
double eps = 1.E-12;
if (IsA() == TH2F::Class() ) eps = 1.E-6;
if (fTsumw != 0 && TMath::Abs( fTsumw - totcont) < TMath::Abs(fTsumw) * eps)
reuseStats = true;
}
if (ncuts) reuseStats = false;
bool reuseEntries = reuseStats;
reuseEntries &= (firstbin==0 && lastbin == inNbin+1);
if (reuseStats) {
Double_t stats[kNstat];
GetStats(stats);
if (!onX) {
stats[2] = stats[4];
stats[3] = stats[5];
}
h1->PutStats(stats);
}
else {
h1->SetEntries( h1->GetEffectiveEntries() );
}
if (reuseEntries) {
h1->SetEntries(fEntries);
}
else {
Double_t entries = TMath::Floor( totcont + 0.5);
if (h1->GetSumw2N()) entries = h1->GetEffectiveEntries();
h1->SetEntries( entries );
}
if (opt.Contains("d")) {
TVirtualPad *padsav = gPad;
TVirtualPad *pad = gROOT->GetSelectedPad();
if (pad) pad->cd();
opt.Remove(opt.First("d"),1);
if (opt.Contains("e")) opt.Remove(opt.First("e"),1);
if (!gPad || !gPad->FindObject(h1)) {
h1->Draw(opt);
} else {
h1->Paint(opt);
}
if (padsav) padsav->cd();
}
return h1;
}
TH1D *TH2::ProjectionX(const char *name, Int_t firstybin, Int_t lastybin, Option_t *option) const
{
return DoProjection(true, name, firstybin, lastybin, option);
}
TH1D *TH2::ProjectionY(const char *name, Int_t firstxbin, Int_t lastxbin, Option_t *option) const
{
return DoProjection(false, name, firstxbin, lastxbin, option);
}
void TH2::PutStats(Double_t *stats)
{
TH1::PutStats(stats);
fTsumwy = stats[4];
fTsumwy2 = stats[5];
fTsumwxy = stats[6];
}
void TH2::Reset(Option_t *option)
{
TH1::Reset(option);
TString opt = option;
opt.ToUpper();
if (opt.Contains("ICE") && !opt.Contains("S")) return;
fTsumwy = 0;
fTsumwy2 = 0;
fTsumwxy = 0;
}
void TH2::SetShowProjectionX(Int_t nbins)
{
GetPainter();
if (fPainter) fPainter->SetShowProjection("x",nbins);
}
void TH2::SetShowProjectionY(Int_t nbins)
{
GetPainter();
if (fPainter) fPainter->SetShowProjection("y",nbins);
}
TH1 *TH2::ShowBackground(Int_t niter, Option_t *option)
{
return (TH1*)gROOT->ProcessLineFast(Form("TSpectrum2::StaticBackground((TH1*)0x%lx,%d,\"%s\")",
(ULong_t)this, niter, option));
}
Int_t TH2::ShowPeaks(Double_t sigma, Option_t *option, Double_t threshold)
{
return (Int_t)gROOT->ProcessLineFast(Form("TSpectrum2::StaticSearch((TH1*)0x%lx,%g,\"%s\",%g)",
(ULong_t)this, sigma, option, threshold));
}
void TH2::Smooth(Int_t ntimes, Option_t *option)
{
Double_t k5a[5][5] = { { 0, 0, 1, 0, 0 },
{ 0, 2, 2, 2, 0 },
{ 1, 2, 5, 2, 1 },
{ 0, 2, 2, 2, 0 },
{ 0, 0, 1, 0, 0 } };
Double_t k5b[5][5] = { { 0, 1, 2, 1, 0 },
{ 1, 2, 4, 2, 1 },
{ 2, 4, 8, 4, 2 },
{ 1, 2, 4, 2, 1 },
{ 0, 1, 2, 1, 0 } };
Double_t k3a[3][3] = { { 0, 1, 0 },
{ 1, 2, 1 },
{ 0, 1, 0 } };
if (ntimes > 1) {
Warning("Smooth","Currently only ntimes=1 is supported");
}
TString opt = option;
opt.ToLower();
Int_t ksize_x=5;
Int_t ksize_y=5;
Double_t *kernel = &k5a[0][0];
if (opt.Contains("k5b")) kernel = &k5b[0][0];
if (opt.Contains("k3a")) {
kernel = &k3a[0][0];
ksize_x=3;
ksize_y=3;
}
Int_t ifirst = fXaxis.GetFirst();
Int_t ilast = fXaxis.GetLast();
Int_t jfirst = fYaxis.GetFirst();
Int_t jlast = fYaxis.GetLast();
Double_t nentries = fEntries;
Int_t nx = GetNbinsX();
Int_t ny = GetNbinsY();
Int_t bufSize = (nx+2)*(ny+2);
Double_t *buf = new Double_t[bufSize];
Double_t *ebuf = 0;
if (fSumw2.fN) ebuf = new Double_t[bufSize];
Int_t i,j,bin;
for (i=ifirst; i<=ilast; i++){
for (j=jfirst; j<=jlast; j++){
bin = GetBin(i,j);
buf[bin] =GetBinContent(bin);
if (ebuf) ebuf[bin]=GetBinError(bin);
}
}
Int_t x_push = (ksize_x-1)/2;
Int_t y_push = (ksize_y-1)/2;
for (i=ifirst; i<=ilast; i++){
for (j=jfirst; j<=jlast; j++) {
Double_t content = 0.0;
Double_t error = 0.0;
Double_t norm = 0.0;
for (Int_t n=0; n<ksize_x; n++) {
for (Int_t m=0; m<ksize_y; m++) {
Int_t xb = i+(n-x_push);
Int_t yb = j+(m-y_push);
if ( (xb >= 1) && (xb <= nx) && (yb >= 1) && (yb <= ny) ) {
bin = GetBin(xb,yb);
Double_t k = kernel[n*ksize_y +m];
if ( k != 0.0 ) {
norm += k;
content += k*buf[bin];
if (ebuf) error += k*k*ebuf[bin]*ebuf[bin];
}
}
}
}
if ( norm != 0.0 ) {
SetBinContent(i,j,content/norm);
if (ebuf) {
error /= (norm*norm);
SetBinError(i,j,sqrt(error));
}
}
}
}
fEntries = nentries;
delete [] buf;
delete [] ebuf;
}
void TH2::Streamer(TBuffer &R__b)
{
if (R__b.IsReading()) {
UInt_t R__s, R__c;
Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
if (R__v > 2) {
R__b.ReadClassBuffer(TH2::Class(), this, R__v, R__s, R__c);
return;
}
TH1::Streamer(R__b);
R__b >> fScalefactor;
R__b >> fTsumwy;
R__b >> fTsumwy2;
R__b >> fTsumwxy;
} else {
R__b.WriteClassBuffer(TH2::Class(),this);
}
}
ClassImp(TH2C)
TH2C::TH2C(): TH2(), TArrayC()
{
SetBinsLength(9);
if (fgDefaultSumw2) Sumw2();
}
TH2C::~TH2C()
{
}
TH2C::TH2C(const char *name,const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
,Int_t nbinsy,Double_t ylow,Double_t yup)
:TH2(name,title,nbinsx,xlow,xup,nbinsy,ylow,yup)
{
TArrayC::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
if (xlow >= xup || ylow >= yup) SetBuffer(fgBufferSize);
}
TH2C::TH2C(const char *name,const char *title,Int_t nbinsx,const Double_t *xbins
,Int_t nbinsy,Double_t ylow,Double_t yup)
:TH2(name,title,nbinsx,xbins,nbinsy,ylow,yup)
{
TArrayC::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
}
TH2C::TH2C(const char *name,const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
,Int_t nbinsy,const Double_t *ybins)
:TH2(name,title,nbinsx,xlow,xup,nbinsy,ybins)
{
TArrayC::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
}
TH2C::TH2C(const char *name,const char *title,Int_t nbinsx,const Double_t *xbins
,Int_t nbinsy,const Double_t *ybins)
:TH2(name,title,nbinsx,xbins,nbinsy,ybins)
{
TArrayC::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
}
TH2C::TH2C(const char *name,const char *title,Int_t nbinsx,const Float_t *xbins
,Int_t nbinsy,const Float_t *ybins)
:TH2(name,title,nbinsx,xbins,nbinsy,ybins)
{
TArrayC::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
}
TH2C::TH2C(const TH2C &h2c) : TH2(), TArrayC()
{
((TH2C&)h2c).Copy(*this);
}
void TH2C::AddBinContent(Int_t bin)
{
if (fArray[bin] < 127) fArray[bin]++;
}
void TH2C::AddBinContent(Int_t bin, Double_t w)
{
Int_t newval = fArray[bin] + Int_t(w);
if (newval > -128 && newval < 128) {fArray[bin] = Char_t(newval); return;}
if (newval < -127) fArray[bin] = -127;
if (newval > 127) fArray[bin] = 127;
}
void TH2C::Copy(TObject &newth2) const
{
TH2::Copy((TH2C&)newth2);
}
TH1 *TH2C::DrawCopy(Option_t *option) const
{
TString opt = option;
opt.ToLower();
if (gPad && !opt.Contains("same")) gPad->Clear();
TH2C *newth2 = (TH2C*)Clone();
newth2->SetDirectory(0);
newth2->SetBit(kCanDelete);
newth2->AppendPad(option);
return newth2;
}
Double_t TH2C::GetBinContent(Int_t bin) const
{
if (fBuffer) ((TH2C*)this)->BufferEmpty();
if (bin < 0) bin = 0;
if (bin >= fNcells) bin = fNcells-1;
if (!fArray) return 0;
return Double_t (fArray[bin]);
}
void TH2C::Reset(Option_t *option)
{
TH2::Reset(option);
TArrayC::Reset();
}
void TH2C::SetBinContent(Int_t bin, Double_t content)
{
fEntries++;
fTsumw = 0;
if (bin < 0) return;
if (bin >= fNcells) return;
fArray[bin] = Char_t (content);
}
void TH2C::SetBinsLength(Int_t n)
{
if (n < 0) n = (fXaxis.GetNbins()+2)*(fYaxis.GetNbins()+2);
fNcells = n;
TArrayC::Set(n);
}
void TH2C::Streamer(TBuffer &R__b)
{
if (R__b.IsReading()) {
UInt_t R__s, R__c;
Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
if (R__v > 2) {
R__b.ReadClassBuffer(TH2C::Class(), this, R__v, R__s, R__c);
return;
}
if (R__v < 2) {
R__b.ReadVersion();
TH1::Streamer(R__b);
TArrayC::Streamer(R__b);
R__b.ReadVersion();
R__b >> fScalefactor;
R__b >> fTsumwy;
R__b >> fTsumwy2;
R__b >> fTsumwxy;
} else {
TH2::Streamer(R__b);
TArrayC::Streamer(R__b);
R__b.CheckByteCount(R__s, R__c, TH2C::IsA());
}
} else {
R__b.WriteClassBuffer(TH2C::Class(),this);
}
}
TH2C& TH2C::operator=(const TH2C &h1)
{
if (this != &h1) ((TH2C&)h1).Copy(*this);
return *this;
}
TH2C operator*(Float_t c1, TH2C &h1)
{
TH2C hnew = h1;
hnew.Scale(c1);
hnew.SetDirectory(0);
return hnew;
}
TH2C operator+(TH2C &h1, TH2C &h2)
{
TH2C hnew = h1;
hnew.Add(&h2,1);
hnew.SetDirectory(0);
return hnew;
}
TH2C operator-(TH2C &h1, TH2C &h2)
{
TH2C hnew = h1;
hnew.Add(&h2,-1);
hnew.SetDirectory(0);
return hnew;
}
TH2C operator*(TH2C &h1, TH2C &h2)
{
TH2C hnew = h1;
hnew.Multiply(&h2);
hnew.SetDirectory(0);
return hnew;
}
TH2C operator/(TH2C &h1, TH2C &h2)
{
TH2C hnew = h1;
hnew.Divide(&h2);
hnew.SetDirectory(0);
return hnew;
}
ClassImp(TH2S)
TH2S::TH2S(): TH2(), TArrayS()
{
SetBinsLength(9);
if (fgDefaultSumw2) Sumw2();
}
TH2S::~TH2S()
{
}
TH2S::TH2S(const char *name,const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
,Int_t nbinsy,Double_t ylow,Double_t yup)
:TH2(name,title,nbinsx,xlow,xup,nbinsy,ylow,yup)
{
TArrayS::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
if (xlow >= xup || ylow >= yup) SetBuffer(fgBufferSize);
}
TH2S::TH2S(const char *name,const char *title,Int_t nbinsx,const Double_t *xbins
,Int_t nbinsy,Double_t ylow,Double_t yup)
:TH2(name,title,nbinsx,xbins,nbinsy,ylow,yup)
{
TArrayS::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
}
TH2S::TH2S(const char *name,const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
,Int_t nbinsy,const Double_t *ybins)
:TH2(name,title,nbinsx,xlow,xup,nbinsy,ybins)
{
TArrayS::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
}
TH2S::TH2S(const char *name,const char *title,Int_t nbinsx,const Double_t *xbins
,Int_t nbinsy,const Double_t *ybins)
:TH2(name,title,nbinsx,xbins,nbinsy,ybins)
{
TArrayS::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
}
TH2S::TH2S(const char *name,const char *title,Int_t nbinsx,const Float_t *xbins
,Int_t nbinsy,const Float_t *ybins)
:TH2(name,title,nbinsx,xbins,nbinsy,ybins)
{
TArrayS::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
}
TH2S::TH2S(const TH2S &h2s) : TH2(), TArrayS()
{
((TH2S&)h2s).Copy(*this);
}
void TH2S::AddBinContent(Int_t bin)
{
if (fArray[bin] < 32767) fArray[bin]++;
}
void TH2S::AddBinContent(Int_t bin, Double_t w)
{
Int_t newval = fArray[bin] + Int_t(w);
if (newval > -32768 && newval < 32768) {fArray[bin] = Short_t(newval); return;}
if (newval < -32767) fArray[bin] = -32767;
if (newval > 32767) fArray[bin] = 32767;
}
void TH2S::Copy(TObject &newth2) const
{
TH2::Copy((TH2S&)newth2);
}
TH1 *TH2S::DrawCopy(Option_t *option) const
{
TString opt = option;
opt.ToLower();
if (gPad && !opt.Contains("same")) gPad->Clear();
TH2S *newth2 = (TH2S*)Clone();
newth2->SetDirectory(0);
newth2->SetBit(kCanDelete);
newth2->AppendPad(option);
return newth2;
}
Double_t TH2S::GetBinContent(Int_t bin) const
{
if (fBuffer) ((TH2C*)this)->BufferEmpty();
if (bin < 0) bin = 0;
if (bin >= fNcells) bin = fNcells-1;
if (!fArray) return 0;
return Double_t (fArray[bin]);
}
void TH2S::Reset(Option_t *option)
{
TH2::Reset(option);
TArrayS::Reset();
}
void TH2S::SetBinContent(Int_t bin, Double_t content)
{
fEntries++;
fTsumw = 0;
if (bin < 0) return;
if (bin >= fNcells) return;
fArray[bin] = Short_t (content);
}
void TH2S::SetBinsLength(Int_t n)
{
if (n < 0) n = (fXaxis.GetNbins()+2)*(fYaxis.GetNbins()+2);
fNcells = n;
TArrayS::Set(n);
}
void TH2S::Streamer(TBuffer &R__b)
{
if (R__b.IsReading()) {
UInt_t R__s, R__c;
Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
if (R__v > 2) {
R__b.ReadClassBuffer(TH2S::Class(), this, R__v, R__s, R__c);
return;
}
if (R__v < 2) {
R__b.ReadVersion();
TH1::Streamer(R__b);
TArrayS::Streamer(R__b);
R__b.ReadVersion();
R__b >> fScalefactor;
R__b >> fTsumwy;
R__b >> fTsumwy2;
R__b >> fTsumwxy;
} else {
TH2::Streamer(R__b);
TArrayS::Streamer(R__b);
R__b.CheckByteCount(R__s, R__c, TH2S::IsA());
}
} else {
R__b.WriteClassBuffer(TH2S::Class(),this);
}
}
TH2S& TH2S::operator=(const TH2S &h1)
{
if (this != &h1) ((TH2S&)h1).Copy(*this);
return *this;
}
TH2S operator*(Float_t c1, TH2S &h1)
{
TH2S hnew = h1;
hnew.Scale(c1);
hnew.SetDirectory(0);
return hnew;
}
TH2S operator+(TH2S &h1, TH2S &h2)
{
TH2S hnew = h1;
hnew.Add(&h2,1);
hnew.SetDirectory(0);
return hnew;
}
TH2S operator-(TH2S &h1, TH2S &h2)
{
TH2S hnew = h1;
hnew.Add(&h2,-1);
hnew.SetDirectory(0);
return hnew;
}
TH2S operator*(TH2S &h1, TH2S &h2)
{
TH2S hnew = h1;
hnew.Multiply(&h2);
hnew.SetDirectory(0);
return hnew;
}
TH2S operator/(TH2S &h1, TH2S &h2)
{
TH2S hnew = h1;
hnew.Divide(&h2);
hnew.SetDirectory(0);
return hnew;
}
ClassImp(TH2I)
TH2I::TH2I(): TH2(), TArrayI()
{
SetBinsLength(9);
if (fgDefaultSumw2) Sumw2();
}
TH2I::~TH2I()
{
}
TH2I::TH2I(const char *name,const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
,Int_t nbinsy,Double_t ylow,Double_t yup)
:TH2(name,title,nbinsx,xlow,xup,nbinsy,ylow,yup)
{
TArrayI::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
if (xlow >= xup || ylow >= yup) SetBuffer(fgBufferSize);
}
TH2I::TH2I(const char *name,const char *title,Int_t nbinsx,const Double_t *xbins
,Int_t nbinsy,Double_t ylow,Double_t yup)
:TH2(name,title,nbinsx,xbins,nbinsy,ylow,yup)
{
TArrayI::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
}
TH2I::TH2I(const char *name,const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
,Int_t nbinsy,const Double_t *ybins)
:TH2(name,title,nbinsx,xlow,xup,nbinsy,ybins)
{
TArrayI::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
}
TH2I::TH2I(const char *name,const char *title,Int_t nbinsx,const Double_t *xbins
,Int_t nbinsy,const Double_t *ybins)
:TH2(name,title,nbinsx,xbins,nbinsy,ybins)
{
TArrayI::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
}
TH2I::TH2I(const char *name,const char *title,Int_t nbinsx,const Float_t *xbins
,Int_t nbinsy,const Float_t *ybins)
:TH2(name,title,nbinsx,xbins,nbinsy,ybins)
{
TArrayI::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
}
TH2I::TH2I(const TH2I &h2i) : TH2(), TArrayI()
{
((TH2I&)h2i).Copy(*this);
}
void TH2I::AddBinContent(Int_t bin)
{
if (fArray[bin] < 2147483647) fArray[bin]++;
}
void TH2I::AddBinContent(Int_t bin, Double_t w)
{
Int_t newval = fArray[bin] + Int_t(w);
if (newval > -2147483647 && newval < 2147483647) {fArray[bin] = Int_t(newval); return;}
if (newval < -2147483647) fArray[bin] = -2147483647;
if (newval > 2147483647) fArray[bin] = 2147483647;
}
void TH2I::Copy(TObject &newth2) const
{
TH2::Copy((TH2I&)newth2);
}
TH1 *TH2I::DrawCopy(Option_t *option) const
{
TString opt = option;
opt.ToLower();
if (gPad && !opt.Contains("same")) gPad->Clear();
TH2I *newth2 = (TH2I*)Clone();
newth2->SetDirectory(0);
newth2->SetBit(kCanDelete);
newth2->AppendPad(option);
return newth2;
}
Double_t TH2I::GetBinContent(Int_t bin) const
{
if (fBuffer) ((TH2C*)this)->BufferEmpty();
if (bin < 0) bin = 0;
if (bin >= fNcells) bin = fNcells-1;
if (!fArray) return 0;
return Double_t (fArray[bin]);
}
void TH2I::Reset(Option_t *option)
{
TH2::Reset(option);
TArrayI::Reset();
}
void TH2I::SetBinContent(Int_t bin, Double_t content)
{
fEntries++;
fTsumw = 0;
if (bin < 0) return;
if (bin >= fNcells) return;
fArray[bin] = Int_t (content);
}
void TH2I::SetBinsLength(Int_t n)
{
if (n < 0) n = (fXaxis.GetNbins()+2)*(fYaxis.GetNbins()+2);
fNcells = n;
TArrayI::Set(n);
}
TH2I& TH2I::operator=(const TH2I &h1)
{
if (this != &h1) ((TH2I&)h1).Copy(*this);
return *this;
}
TH2I operator*(Float_t c1, TH2I &h1)
{
TH2I hnew = h1;
hnew.Scale(c1);
hnew.SetDirectory(0);
return hnew;
}
TH2I operator+(TH2I &h1, TH2I &h2)
{
TH2I hnew = h1;
hnew.Add(&h2,1);
hnew.SetDirectory(0);
return hnew;
}
TH2I operator-(TH2I &h1, TH2I &h2)
{
TH2I hnew = h1;
hnew.Add(&h2,-1);
hnew.SetDirectory(0);
return hnew;
}
TH2I operator*(TH2I &h1, TH2I &h2)
{
TH2I hnew = h1;
hnew.Multiply(&h2);
hnew.SetDirectory(0);
return hnew;
}
TH2I operator/(TH2I &h1, TH2I &h2)
{
TH2I hnew = h1;
hnew.Divide(&h2);
hnew.SetDirectory(0);
return hnew;
}
ClassImp(TH2F)
TH2F::TH2F(): TH2(), TArrayF()
{
SetBinsLength(9);
if (fgDefaultSumw2) Sumw2();
}
TH2F::~TH2F()
{
}
TH2F::TH2F(const char *name,const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
,Int_t nbinsy,Double_t ylow,Double_t yup)
:TH2(name,title,nbinsx,xlow,xup,nbinsy,ylow,yup)
{
TArrayF::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
if (xlow >= xup || ylow >= yup) SetBuffer(fgBufferSize);
}
TH2F::TH2F(const char *name,const char *title,Int_t nbinsx,const Double_t *xbins
,Int_t nbinsy,Double_t ylow,Double_t yup)
:TH2(name,title,nbinsx,xbins,nbinsy,ylow,yup)
{
TArrayF::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
}
TH2F::TH2F(const char *name,const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
,Int_t nbinsy,const Double_t *ybins)
:TH2(name,title,nbinsx,xlow,xup,nbinsy,ybins)
{
TArrayF::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
}
TH2F::TH2F(const char *name,const char *title,Int_t nbinsx,const Double_t *xbins
,Int_t nbinsy,const Double_t *ybins)
:TH2(name,title,nbinsx,xbins,nbinsy,ybins)
{
TArrayF::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
}
TH2F::TH2F(const char *name,const char *title,Int_t nbinsx,const Float_t *xbins
,Int_t nbinsy,const Float_t *ybins)
:TH2(name,title,nbinsx,xbins,nbinsy,ybins)
{
TArrayF::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
}
TH2F::TH2F(const TMatrixFBase &m)
:TH2("TMatrixFBase","",m.GetNcols(),m.GetColLwb(),1+m.GetColUpb(),m.GetNrows(),m.GetRowLwb(),1+m.GetRowUpb())
{
TArrayF::Set(fNcells);
Int_t ilow = m.GetRowLwb();
Int_t iup = m.GetRowUpb();
Int_t jlow = m.GetColLwb();
Int_t jup = m.GetColUpb();
for (Int_t i=ilow;i<=iup;i++) {
for (Int_t j=jlow;j<=jup;j++) {
SetCellContent(j-jlow+1,i-ilow+1,m(i,j));
}
}
}
TH2F::TH2F(const TH2F &h2f) : TH2(), TArrayF()
{
((TH2F&)h2f).Copy(*this);
}
void TH2F::Copy(TObject &newth2) const
{
TH2::Copy((TH2F&)newth2);
}
TH1 *TH2F::DrawCopy(Option_t *option) const
{
TString opt = option;
opt.ToLower();
if (gPad && !opt.Contains("same")) gPad->Clear();
TH2F *newth2 = (TH2F*)Clone();
newth2->SetDirectory(0);
newth2->SetBit(kCanDelete);
newth2->AppendPad(option);
return newth2;
}
Double_t TH2F::GetBinContent(Int_t bin) const
{
if (fBuffer) ((TH2C*)this)->BufferEmpty();
if (bin < 0) bin = 0;
if (bin >= fNcells) bin = fNcells-1;
if (!fArray) return 0;
return Double_t (fArray[bin]);
}
void TH2F::Reset(Option_t *option)
{
TH2::Reset(option);
TArrayF::Reset();
}
void TH2F::SetBinContent(Int_t bin, Double_t content)
{
fEntries++;
fTsumw = 0;
if (bin < 0) return;
if (bin >= fNcells) return;
fArray[bin] = Float_t (content);
}
void TH2F::SetBinsLength(Int_t n)
{
if (n < 0) n = (fXaxis.GetNbins()+2)*(fYaxis.GetNbins()+2);
fNcells = n;
TArrayF::Set(n);
}
void TH2F::Streamer(TBuffer &R__b)
{
if (R__b.IsReading()) {
UInt_t R__s, R__c;
Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
if (R__v > 2) {
R__b.ReadClassBuffer(TH2F::Class(), this, R__v, R__s, R__c);
return;
}
if (R__v < 2) {
R__b.ReadVersion();
TH1::Streamer(R__b);
TArrayF::Streamer(R__b);
R__b.ReadVersion();
R__b >> fScalefactor;
R__b >> fTsumwy;
R__b >> fTsumwy2;
R__b >> fTsumwxy;
} else {
TH2::Streamer(R__b);
TArrayF::Streamer(R__b);
R__b.CheckByteCount(R__s, R__c, TH2F::IsA());
}
} else {
R__b.WriteClassBuffer(TH2F::Class(),this);
}
}
TH2F& TH2F::operator=(const TH2F &h1)
{
if (this != &h1) ((TH2F&)h1).Copy(*this);
return *this;
}
TH2F operator*(Float_t c1, TH2F &h1)
{
TH2F hnew = h1;
hnew.Scale(c1);
hnew.SetDirectory(0);
return hnew;
}
TH2F operator*(TH2F &h1, Float_t c1)
{
TH2F hnew = h1;
hnew.Scale(c1);
hnew.SetDirectory(0);
return hnew;
}
TH2F operator+(TH2F &h1, TH2F &h2)
{
TH2F hnew = h1;
hnew.Add(&h2,1);
hnew.SetDirectory(0);
return hnew;
}
TH2F operator-(TH2F &h1, TH2F &h2)
{
TH2F hnew = h1;
hnew.Add(&h2,-1);
hnew.SetDirectory(0);
return hnew;
}
TH2F operator*(TH2F &h1, TH2F &h2)
{
TH2F hnew = h1;
hnew.Multiply(&h2);
hnew.SetDirectory(0);
return hnew;
}
TH2F operator/(TH2F &h1, TH2F &h2)
{
TH2F hnew = h1;
hnew.Divide(&h2);
hnew.SetDirectory(0);
return hnew;
}
ClassImp(TH2D)
TH2D::TH2D(): TH2(), TArrayD()
{
SetBinsLength(9);
if (fgDefaultSumw2) Sumw2();
}
TH2D::~TH2D()
{
}
TH2D::TH2D(const char *name,const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
,Int_t nbinsy,Double_t ylow,Double_t yup)
:TH2(name,title,nbinsx,xlow,xup,nbinsy,ylow,yup)
{
TArrayD::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
if (xlow >= xup || ylow >= yup) SetBuffer(fgBufferSize);
}
TH2D::TH2D(const char *name,const char *title,Int_t nbinsx,const Double_t *xbins
,Int_t nbinsy,Double_t ylow,Double_t yup)
:TH2(name,title,nbinsx,xbins,nbinsy,ylow,yup)
{
TArrayD::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
}
TH2D::TH2D(const char *name,const char *title,Int_t nbinsx,Double_t xlow,Double_t xup
,Int_t nbinsy,const Double_t *ybins)
:TH2(name,title,nbinsx,xlow,xup,nbinsy,ybins)
{
TArrayD::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
}
TH2D::TH2D(const char *name,const char *title,Int_t nbinsx,const Double_t *xbins
,Int_t nbinsy,const Double_t *ybins)
:TH2(name,title,nbinsx,xbins,nbinsy,ybins)
{
TArrayD::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
}
TH2D::TH2D(const char *name,const char *title,Int_t nbinsx,const Float_t *xbins
,Int_t nbinsy,const Float_t *ybins)
:TH2(name,title,nbinsx,xbins,nbinsy,ybins)
{
TArrayD::Set(fNcells);
if (fgDefaultSumw2) Sumw2();
}
TH2D::TH2D(const TMatrixDBase &m)
:TH2("TMatrixDBase","",m.GetNcols(),m.GetColLwb(),1+m.GetColUpb(),m.GetNrows(),m.GetRowLwb(),1+m.GetRowUpb())
{
TArrayD::Set(fNcells);
Int_t ilow = m.GetRowLwb();
Int_t iup = m.GetRowUpb();
Int_t jlow = m.GetColLwb();
Int_t jup = m.GetColUpb();
for (Int_t i=ilow;i<=iup;i++) {
for (Int_t j=jlow;j<=jup;j++) {
SetCellContent(j-jlow+1,i-ilow+1,m(i,j));
}
}
if (fgDefaultSumw2) Sumw2();
}
TH2D::TH2D(const TH2D &h2d) : TH2(), TArrayD()
{
((TH2D&)h2d).Copy(*this);
}
void TH2D::Copy(TObject &newth2) const
{
TH2::Copy((TH2D&)newth2);
}
TH1 *TH2D::DrawCopy(Option_t *option) const
{
TString opt = option;
opt.ToLower();
if (gPad && !opt.Contains("same")) gPad->Clear();
TH2D *newth2 = (TH2D*)Clone();
newth2->SetDirectory(0);
newth2->SetBit(kCanDelete);
newth2->AppendPad(option);
return newth2;
}
Double_t TH2D::GetBinContent(Int_t bin) const
{
if (fBuffer) ((TH2C*)this)->BufferEmpty();
if (bin < 0) bin = 0;
if (bin >= fNcells) bin = fNcells-1;
if (!fArray) return 0;
return Double_t (fArray[bin]);
}
void TH2D::Reset(Option_t *option)
{
TH2::Reset(option);
TArrayD::Reset();
}
void TH2D::SetBinContent(Int_t bin, Double_t content)
{
fEntries++;
fTsumw = 0;
if (bin < 0) return;
if (bin >= fNcells) return;
fArray[bin] = Double_t (content);
}
void TH2D::SetBinsLength(Int_t n)
{
if (n < 0) n = (fXaxis.GetNbins()+2)*(fYaxis.GetNbins()+2);
fNcells = n;
TArrayD::Set(n);
}
void TH2D::Streamer(TBuffer &R__b)
{
if (R__b.IsReading()) {
UInt_t R__s, R__c;
Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
if (R__v > 2) {
R__b.ReadClassBuffer(TH2D::Class(), this, R__v, R__s, R__c);
return;
}
if (R__v < 2) {
R__b.ReadVersion();
TH1::Streamer(R__b);
TArrayD::Streamer(R__b);
R__b.ReadVersion();
R__b >> fScalefactor;
R__b >> fTsumwy;
R__b >> fTsumwy2;
R__b >> fTsumwxy;
} else {
TH2::Streamer(R__b);
TArrayD::Streamer(R__b);
R__b.CheckByteCount(R__s, R__c, TH2D::IsA());
}
} else {
R__b.WriteClassBuffer(TH2D::Class(),this);
}
}
TH2D& TH2D::operator=(const TH2D &h1)
{
if (this != &h1) ((TH2D&)h1).Copy(*this);
return *this;
}
TH2D operator*(Float_t c1, TH2D &h1)
{
TH2D hnew = h1;
hnew.Scale(c1);
hnew.SetDirectory(0);
return hnew;
}
TH2D operator+(TH2D &h1, TH2D &h2)
{
TH2D hnew = h1;
hnew.Add(&h2,1);
hnew.SetDirectory(0);
return hnew;
}
TH2D operator-(TH2D &h1, TH2D &h2)
{
TH2D hnew = h1;
hnew.Add(&h2,-1);
hnew.SetDirectory(0);
return hnew;
}
TH2D operator*(TH2D &h1, TH2D &h2)
{
TH2D hnew = h1;
hnew.Multiply(&h2);
hnew.SetDirectory(0);
return hnew;
}
TH2D operator/(TH2D &h1, TH2D &h2)
{
TH2D hnew = h1;
hnew.Divide(&h2);
hnew.SetDirectory(0);
return hnew;
}