/************************************************************************* * Copyright(c) 1998-2048, ALICE Experiment at CERN, All rights reserved. * * * * Author: The ALICE Off-line Project. * * Contributors are mentioned in the code where appropriate. * * * * Permission to use, copy, modify and distribute this software and its * * documentation strictly for non-commercial purposes is hereby granted * * without fee, provided that the above copyright notice appears in all * * copies and that both the copyright notice and this permission notice * * appear in the supporting documentation. The authors make no claims * * about the suitability of this software for any purpose. It is * * provided "as is" without express or implied warranty. * **************************************************************************/ // Author: Dariusz Miskowiec 2007 //============================================================================= // multidimensional histogram // Physically, the data are kept in one single one-dimensional histogram. // Projecting on 1, 2, and n-1 dimensions is implemented. // The histogram can be saved on root file as the one-dimensional data // histogram and the axes, thus eternal forward compatibility is ensured. //============================================================================= #include #include #include #include #include #include #include "AliUnicorHN.h" ClassImp(AliUnicorHN) //============================================================================= AliUnicorHN::AliUnicorHN(Char_t *nam, Int_t ndim, TAxis **ax) : TH1D(nam, nam, Albins(ndim,ax), 0.5, Albins(ndim,ax)+0.5), fNdim(ndim) { // Constructor for building from scratch. // Above, we just have managed to create a single one-dimensional histogram // with number of bins equal to the product of the numbers of bins in all // dimensions. For easy inspection the histogram range was set to -0.5,n-0.5. for (int i=0; iSetName(Form("axis%d",i)); for (int i=0; iCopy(fAxis[i]); // for speed reasons, number of bins of each axis is stored on fNbins // and the running product of the latter on fMbins // so fMbins = {...,fNbins[fNdim-2]*fNbins[fNdim-1],fNbins[fNdim-1],1} fMbins[fNdim-1] = 1; for (int i=0; i0; i--) fMbins[i-1] = fMbins[i]*fNbins[i]; printf(" %d-dimensional histogram %s with %d bins created\n",fNdim,nam,GetNbinsX()); } //============================================================================= AliUnicorHN::AliUnicorHN(Char_t *filnam, Char_t *nam) : TH1D(*((TH1D*) TFile::Open(filnam,"read")->GetDirectory(nam)->Get("histo"))), fNdim(0) { // Constructor for reading from file. TFile *f = TFile::Open(filnam,"read"); f->cd(nam); TAxis *ax[fgkMaxNdim]; for (fNdim=0; fNdimGet(Form("axis%d",fNdim)); if (ax[fNdim]) ax[fNdim]->Copy(fAxis[fNdim]); else break; } f->Close(); fMbins[fNdim-1] = 1; for (int i=0; i0; i--) fMbins[i-1] = fMbins[i]*fNbins[i]; if (GetNbinsX()!=Albins(fNdim,ax)) { printf("number of bins of histo %d differs from product of nbins of axes %d\n", GetNbinsX(),Albins(fNdim,ax)); printf("bombing\n"); exit(-1); } printf("%d-dimensional histogram %s with %d bins read from file %s\n", fNdim,nam,GetNbinsX(),filnam); } //============================================================================= Int_t AliUnicorHN::Albins(Int_t n, TAxis **ax) { // Calculate product of nbins of ax[0]...ax[n-1] // (= total number of bins of the multidimensional histogram to be made). Int_t nbins = 1; // while (n--) nbins *= ax[n]->GetNbins(); for (int i=0; iGetNbins(); return nbins; } //============================================================================= Int_t AliUnicorHN::MulToOne(const Int_t * const k) const { // Calculate the 1-dim index n from n-dim indices k[fNdim]. // Valid k[i] should be between 0 and fNbins[i]-1. // Valid result will be between 0 and GetNbinsX()-1. // Return -1 if under- or overflow in any dimension. Int_t n = 0; for (int i=0; i=fNbins[i]) return -1; n += fMbins[i]*k[i]; } return n; } //============================================================================= Int_t AliUnicorHN::MulToOne(Double_t *x) { // Calculate the 1-dim index n from n-dim vector x, representing the // abscissa of the n-dim histogram. The result will be between 0 and // GetNbinsX()-1. Int_t k[fgkMaxNdim]; for (int i=0; imkdir(GetName()); if (dest) { dest->cd(); nbytes += histo.Write("histo"); for (int i=0; iGetPath()); dest->cd(".."); } return nbytes; } //============================================================================= AliUnicorHN *AliUnicorHN::ProjectAlong(char *nam, Int_t dim, Int_t first, Int_t last) { // Reduce dimension dim by summing up its bins between first and last. // Use root convention: bin=1 is the first bin, bin=nbins is the last. // Return the resulting fNdim-1 dimensional histogram. if (dim<0 || dim>fNdim-1) return 0; if (first<0) first = 1; if (last<0) last = fNbins[dim]; // create new (reduced) histogram TAxis *ax[fgkMaxNdim-1]; int n=0; for (int i=0; ilast) continue; n = 0; for (int j=0; jMulToOne(m); his->AddBinContent(n+1,GetBinContent(i+1)); eis->AddBinContent(n+1,GetBinError(i+1)*GetBinError(i+1)); } // combine content and errors in one histogram for (int i=0; iGetNbinsX(); i++) { his->SetBinError(i+1,sqrt(eis->GetBinContent(i+1))); } his->SetLineColor(this->GetLineColor()); his->SetFillColor(this->GetFillColor()); his->SetMarkerColor(this->GetMarkerColor()); his->SetMarkerStyle(this->GetMarkerStyle()); // some cleanup delete eis; return his; } //============================================================================= TH1D *AliUnicorHN::ProjectOn(char *nam, Int_t dim, const Int_t * const first, const Int_t * const last) const { // Project on dimension dim. Use only bins between first[i] and last[i]. // Use root convention: bin=1 is the first bin, bin=nbins is the last. // first[i]=0 means from the first bin // last[i]=0 means till the last bin if (dim<0 || dim>fNdim-1) return 0; // calculate the projection; lowest index 0 double *yy = new double[fNbins[dim]]; // value double *ey = new double[fNbins[dim]]; // error for (int i=0; ilast[j]) {isgood = 0; break;} } if (isgood) { yy[k[dim]]+=GetBinContent(i+1); ey[k[dim]]+=GetBinError(i+1)*GetBinError(i+1); } } // make the projection histogram // use name nam if specified; otherwise generate one TH1D *his; char *name = strlen(nam)? nam : Form("%s_proj%d",GetName(),dim); if (fAxis[dim].IsVariableBinSize()) his = new TH1D(name,name,fNbins[dim],fAxis[dim].GetXbins()->GetArray()); else his = new TH1D(name,name,fNbins[dim],fAxis[dim].GetXmin(),fAxis[dim].GetXmax()); his->SetXTitle(fAxis[dim].GetTitle()); // or his->GetXaxis()->ImportAttributes(ax); his->Sumw2(); his->SetLineColor(this->GetLineColor()); his->SetFillColor(this->GetFillColor()); his->SetMarkerColor(this->GetMarkerColor()); his->SetMarkerStyle(this->GetMarkerStyle()); for (int i=0; iGetNbinsX(); i++) { his->SetBinContent(i+1,yy[i]); his->SetBinError(i+1,sqrt(ey[i])); } // some cleanup delete [] yy; delete [] ey; delete [] k; // if (name!=nam) delete [] name; return his; } //============================================================================= TH1D *AliUnicorHN::ProjectOn(char *nam, Int_t dim, const Double_t * const first, const Double_t * const last) { // Project on dimension dim. Use only bins between first[i] and last[i]. if (dim<0 || dim>fNdim-1) return 0; Int_t kfirst[fgkMaxNdim]; Int_t klast[fgkMaxNdim]; for (int i=0; ifNdim-1) return 0; if (dim1<0 || dim1>fNdim-1) return 0; // calculate the projection double **yy = new double*[fNbins[dim0]]; // value double **ey = new double*[fNbins[dim0]]; // error for (int i=0; ilast[j]) {isgood = 0; break;} } if (isgood) { yy[k[dim0]][k[dim1]]+=GetBinContent(i+1); ey[k[dim0]][k[dim1]]+=GetBinError(i+1)*GetBinError(i+1); } } // make the projection histogram // use name nam if specified; otherwise generate one TH2D *his=0; char *name = strlen(nam)? nam : Form("%s_proj%dvs%d",GetName(),dim1,dim0); if (fAxis[dim0].IsVariableBinSize() && fAxis[dim1].IsVariableBinSize()) his = new TH2D(name,name, fNbins[dim0],fAxis[dim0].GetXbins()->GetArray(), fNbins[dim1],fAxis[dim1].GetXbins()->GetArray()); else if (!fAxis[dim0].IsVariableBinSize() && fAxis[dim1].IsVariableBinSize()) his = new TH2D(name,name, fNbins[dim0],fAxis[dim0].GetXmin(),fAxis[dim0].GetXmax(), fNbins[dim1],fAxis[dim1].GetXbins()->GetArray()); else if (fAxis[dim0].IsVariableBinSize() && !fAxis[dim1].IsVariableBinSize()) his = new TH2D(name,name, fNbins[dim0],fAxis[dim0].GetXbins()->GetArray(), fNbins[dim1],fAxis[dim1].GetXmin(),fAxis[dim1].GetXmax()); else if (!fAxis[dim0].IsVariableBinSize() && !fAxis[dim1].IsVariableBinSize()) his = new TH2D(name,name, fNbins[dim0],fAxis[dim0].GetXmin(),fAxis[dim0].GetXmax(), fNbins[dim1],fAxis[dim1].GetXmin(),fAxis[dim1].GetXmax()); his->SetXTitle(fAxis[dim0].GetTitle()); his->SetYTitle(fAxis[dim1].GetTitle()); // or his->GetXaxis()->ImportAttributes(ax0); // or his->GetYaxis()->ImportAttributes(ax1); his->Sumw2(); his->SetLineColor(this->GetLineColor()); his->SetFillColor(this->GetFillColor()); his->SetMarkerColor(this->GetMarkerColor()); his->SetMarkerStyle(this->GetMarkerStyle()); for (int i=0; iGetNbinsX(); i++) for (int j=0; jGetNbinsY(); j++) { his->SetBinContent(i+1,j+1,yy[i][j]); his->SetBinError(i+1,j+1,sqrt(ey[i][j])); } // some cleanup for (int i=0; i