3 //--------------------------------------------------------------------//
5 // AliCFEffGrid Class //
6 // Class to handle efficiency grids //
8 // -- Author : S.Arcelli //
12 //--------------------------------------------------------------------//
19 #include "AliCFEffGrid.h"
21 //____________________________________________________________________
22 ClassImp(AliCFEffGrid)
24 //____________________________________________________________________
25 AliCFEffGrid::AliCFEffGrid() :
32 // default constructor
36 //____________________________________________________________________
37 AliCFEffGrid::AliCFEffGrid(const Char_t* name, const Char_t* title, const Int_t nVarIn, const Int_t * nBinIn, const Float_t *binLimitsIn) :
38 AliCFGrid(name,title,nVarIn,nBinIn,binLimitsIn),
48 //____________________________________________________________________
49 AliCFEffGrid::AliCFEffGrid(const Char_t* name, const Char_t* title, const AliCFContainer &c) :
50 AliCFGrid(name,title,c.GetNVar(),c.GetNBins(),c.GetBinLimits()),
59 //assign the container;
62 //____________________________________________________________________
63 AliCFEffGrid::AliCFEffGrid(const AliCFEffGrid& eff) : AliCFGrid(),
71 ((AliCFEffGrid &)eff).Copy(*this);
74 //____________________________________________________________________
75 AliCFEffGrid::~AliCFEffGrid()
82 //____________________________________________________________________
83 AliCFEffGrid &AliCFEffGrid::operator=(const AliCFEffGrid &eff)
89 ((AliCFEffGrid &) eff).Copy(*this);
92 //____________________________________________________________________
94 void AliCFEffGrid::CalculateEfficiency(Int_t istep1,Int_t istep2)
97 // Calculate the efficiency matrix and its error between selection
98 // Steps istep1 and istep2
103 AliCFGrid *num=fContainer->GetGrid(fSelNum);
104 AliCFGrid *den=fContainer->GetGrid(fSelDen);
107 this->Divide(num,den,1.,1.,"B");
109 Int_t nEmptyBinsNum=0;
110 Int_t nEmptyBinsNumAndDen=0;
111 for(Int_t iel=0;iel<fNDim;iel++){
112 if(den->GetElement(iel)>0){
113 if(num->GetElement(iel)==0)nEmptyBinsNum++; //num==0,den!=0
116 nEmptyBinsNumAndDen++;
119 // Some monitoring printout:
120 AliInfo(Form("Efficiency calculated for steps %i and %i: %i empty bins in the numerator && !denominator and %i empty bins in numerator && denominator were found.",fSelNum,fSelDen,nEmptyBinsNumAndDen,nEmptyBinsNum));
121 AliInfo(Form("The correction map contains %i empty bins ",nEmptyBinsNum+nEmptyBinsNumAndDen));
123 //_____________________________________________________________________
124 Float_t AliCFEffGrid::GetAverage() const
127 // Get the average efficiency
133 for(Int_t i=0;i<fNDim;i++){
134 valnum+=fContainer->GetGrid(fSelNum)->GetElement(i);
135 valden+=fContainer->GetGrid(fSelDen)->GetElement(i);
137 if(valden>0)val=valnum/valden;
138 AliInfo(Form(" The Average Efficiency = %f ",val));
142 //_____________________________________________________________________
143 Float_t AliCFEffGrid::GetAverage(Float_t *varMin, Float_t* varMax ) const
146 // Get ave efficiency in a range
151 Int_t *indexMin = new Int_t[fNVar];
152 Int_t *indexMax = new Int_t[fNVar];
153 Int_t *index = new Int_t[fNVar];
155 //Find out the min and max bins
157 for(Int_t i=0;i<fNVar;i++){
158 Float_t xmin=varMin[i]; // the min values
159 Float_t xmax=varMax[i]; // the max values
160 Int_t nbins=fNVarBins[i]+1;
161 Float_t *bins=new Float_t[nbins];
162 for(Int_t ibin =0;ibin<nbins;ibin++){
163 bins[ibin] = fVarBinLimits[ibin+fOffset[i]];
165 indexMin[i] = TMath::BinarySearch(nbins,bins,xmin);
166 indexMax[i] = TMath::BinarySearch(nbins,bins,xmax);
167 if(xmax>=bins[nbins-1]){
168 indexMax[i]=indexMax[i]-1;
175 for(Int_t i=0;i<fNDim;i++){
176 for (Int_t j=0;j<fNVar;j++)index[j]=GetBinIndex(j,i);
178 for (Int_t j=0;j<fNVar;j++){
179 if(!(index[j]>=indexMin[j] && index[j]<=indexMax[j]))isIn=kFALSE;
182 valnum+=fContainer->GetGrid(fSelNum)->GetElement(i);
183 valden+=fContainer->GetGrid(fSelDen)->GetElement(i);
189 if(valden>0)val=valnum/valden;
190 AliInfo(Form(" the Average Efficiency = %f ",val));
193 //____________________________________________________________________
194 void AliCFEffGrid::Copy(TObject& eff) const
200 AliCFEffGrid& target = (AliCFEffGrid &) eff;
202 target.fSelNum=fSelNum;
203 target.fSelDen=fSelDen;
205 target.fContainer=fContainer;
207 //___________________________________________________________________
208 TH1F *AliCFEffGrid::Project(Int_t ivar) const
211 // Make a 1D projection along variable ivar
215 Int_t nbins =fNVarBins[ivar];
216 Float_t *bins = new Float_t[nbins+1];
217 for(Int_t ibin =0;ibin<nbins+1;ibin++){
218 bins[ibin] = fVarBinLimits[ibin+fOffset[ivar]];
222 sprintf(pname,"%s%s%i%i%s%i",GetName(),"_SelStep",fSelNum,fSelDen,"_proj1D_var", ivar);
224 sprintf(htitle,"%s%s%i%i%s%i",GetName(),"_SelStep",fSelNum,fSelDen,"_proj1D_var", ivar);
227 proj1D =new TH1F(pname,htitle, nbins, bins);
230 proj1D->Divide(fContainer->GetGrid(fSelNum)->Project(ivar),fContainer->GetGrid(fSelDen)->Project(ivar),1.,1.,"B");
235 //___________________________________________________________________
236 TH2F *AliCFEffGrid::Project(Int_t ivar1,Int_t ivar2) const
239 // Make a 2D projection along variable ivar1,ivar2
244 Int_t nbins1 =fNVarBins[ivar1];
245 Float_t *bins1 = new Float_t[nbins1+1];
246 Int_t nbins2 =fNVarBins[ivar2];
247 Float_t *bins2 = new Float_t[nbins2+1];
248 for(Int_t ibin1 =0;ibin1<nbins1+1;ibin1++){
249 bins1[ibin1] = fVarBinLimits[ibin1+fOffset[ivar1]];
251 for(Int_t ibin2 =0;ibin2<nbins2+1;ibin2++){
252 bins2[ibin2] = fVarBinLimits[ibin2+fOffset[ivar2]];
256 sprintf(pname,"%s%s%i%i%s%i%i",GetName(),"_SelStep",fSelNum,fSelDen,"_proj2D_var", ivar1,ivar2);
258 sprintf(htitle,"%s%s%i%i%s%i%i",GetName(),"_SelStep",fSelNum,fSelDen,"_proj2D_var",ivar1,ivar2);
261 proj2D =new TH2F(pname,htitle, nbins1,bins1,nbins2,bins2);
264 proj2D->Divide(fContainer->GetGrid(fSelNum)->Project(ivar1,ivar2),fContainer->GetGrid(fSelDen)->Project(ivar1,ivar2),1.,1.,"B");
270 //___________________________________________________________________
271 TH3F *AliCFEffGrid::Project(Int_t ivar1, Int_t ivar2, Int_t ivar3) const
274 // Make a 3D projection along variable ivar1,ivar2,ivar3
279 Int_t nbins1 =fNVarBins[ivar1];
280 Int_t nbins2 =fNVarBins[ivar2];
281 Int_t nbins3 =fNVarBins[ivar3];
283 Float_t *bins1 = new Float_t[nbins1+1];
284 Float_t *bins2 = new Float_t[nbins2+1];
285 Float_t *bins3 = new Float_t[nbins3+1];
287 for(Int_t ibin =0;ibin<nbins1+1;ibin++){
288 bins1[ibin] = fVarBinLimits[ibin+fOffset[ivar1]];
290 for(Int_t ibin =0;ibin<nbins2+1;ibin++){
291 bins2[ibin] = fVarBinLimits[ibin+fOffset[ivar2]];
293 for(Int_t ibin =0;ibin<nbins3+1;ibin++){
294 bins3[ibin] = fVarBinLimits[ibin+fOffset[ivar3]];
298 sprintf(pname,"%s%s%i%i%s%i%i%i",GetName(),"_SelStep",fSelNum,fSelDen,"_proj3D_var",ivar1,ivar2,ivar3);
300 sprintf(htitle,"%s%s%i%i%s%i%i%i",GetName(),"_SelStep",fSelNum,fSelDen,"_proj3D_var",ivar1,ivar2,ivar3);
304 proj3D =new TH3F(pname,htitle, nbins1, bins1,nbins2,bins2,nbins3,bins3);
307 proj3D->Divide(fContainer->GetGrid(fSelNum)->Project(ivar1,ivar2,ivar3),fContainer->GetGrid(fSelDen)->Project(ivar1,ivar2,ivar3),1.,1.,"B");