]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/AliCFEffGrid.cxx
test svn commit, removing
[u/mrichter/AliRoot.git] / CORRFW / AliCFEffGrid.cxx
1 /* $Id$ */
2
3 //--------------------------------------------------------------------//
4 //                                                                    //
5 // AliCFEffGrid Class                                              //
6 // Class to handle efficiency grids                                   // 
7 //                                                                    //
8 // -- Author : S.Arcelli                                              //
9 //                                                                    //
10 //                                                                    //
11 //                                                                    //
12 //--------------------------------------------------------------------//
13 //
14 //
15 #include <TROOT.h>
16 #include <TMath.h>
17 #include <TFile.h>
18 #include <AliLog.h>
19 #include "AliCFEffGrid.h"
20
21 //____________________________________________________________________
22 ClassImp(AliCFEffGrid)
23
24 //____________________________________________________________________
25 AliCFEffGrid::AliCFEffGrid() : 
26   AliCFGrid(),
27   fContainer(0x0),
28   fSelNum(-1),
29   fSelDen(-1)
30 {
31   //
32   // default constructor
33   //
34 }
35
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),
39   fContainer(0x0),
40   fSelNum(-1),
41   fSelDen(-1)
42 {
43   //
44   // ctor
45   //
46   SumW2();
47 }
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()),
51   fContainer(0x0),
52   fSelNum(-1),
53   fSelDen(-1)
54 {
55   //
56   // main constructor
57   //
58   SumW2();
59   //assign the container;
60   fContainer=&c;
61 }
62 //____________________________________________________________________
63 AliCFEffGrid::AliCFEffGrid(const AliCFEffGrid& eff) :   AliCFGrid(),
64   fContainer(0x0),
65   fSelNum(-1),
66   fSelDen(-1)
67 {
68   //
69   // copy constructor
70   //
71   ((AliCFEffGrid &)eff).Copy(*this);
72 }
73
74 //____________________________________________________________________
75 AliCFEffGrid::~AliCFEffGrid()
76 {
77   //
78   // destructor
79   //
80 }
81
82 //____________________________________________________________________
83 AliCFEffGrid &AliCFEffGrid::operator=(const AliCFEffGrid &eff)
84 {
85   //
86   // assigment operator
87   //
88   if (this != &eff)
89     ((AliCFEffGrid &) eff).Copy(*this);
90   return *this;
91
92 //____________________________________________________________________
93
94 void AliCFEffGrid::CalculateEfficiency(Int_t istep1,Int_t istep2)
95 {
96   //
97   // Calculate the efficiency matrix and its error between selection
98   // Steps istep1 and istep2
99   //
100
101   fSelNum=istep1;
102   fSelDen=istep2;
103   AliCFGrid *num=fContainer->GetGrid(fSelNum);
104   AliCFGrid *den=fContainer->GetGrid(fSelDen);
105   num->SumW2();
106   den->SumW2();
107   this->Divide(num,den,1.,1.,"B");
108
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
114       }
115     else{
116       nEmptyBinsNumAndDen++;
117     }
118   }    
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));
122
123 //_____________________________________________________________________
124 Float_t AliCFEffGrid::GetAverage() const 
125 {
126   //
127   // Get the average efficiency 
128   //
129
130   Float_t val=0;
131   Float_t valnum=0;
132   Float_t valden=0;
133   for(Int_t i=0;i<fNDim;i++){
134     valnum+=fContainer->GetGrid(fSelNum)->GetElement(i);
135     valden+=fContainer->GetGrid(fSelDen)->GetElement(i);
136   }
137   if(valden>0)val=valnum/valden;
138   AliInfo(Form(" The Average Efficiency = %f ",val)); 
139
140   return val;
141
142 //_____________________________________________________________________
143 Float_t AliCFEffGrid::GetAverage(Float_t *varMin, Float_t* varMax ) const 
144 {
145   //
146   // Get ave efficiency in a range
147   //
148
149
150   Float_t val=0;
151   Int_t *indexMin = new Int_t[fNVar];
152   Int_t *indexMax = new Int_t[fNVar];
153   Int_t *index    = new Int_t[fNVar];
154   
155   //Find out the min and max bins
156   
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]];
164     }
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;
169     }  
170     delete [] bins;
171   }
172   
173   Float_t valnum=0;
174   Float_t valden=0;
175   for(Int_t i=0;i<fNDim;i++){
176     for (Int_t j=0;j<fNVar;j++)index[j]=GetBinIndex(j,i);
177     Bool_t isIn=kTRUE;
178     for (Int_t j=0;j<fNVar;j++){
179       if(!(index[j]>=indexMin[j] && index[j]<=indexMax[j]))isIn=kFALSE;   
180     }
181     if(isIn){
182       valnum+=fContainer->GetGrid(fSelNum)->GetElement(i);
183       valden+=fContainer->GetGrid(fSelDen)->GetElement(i);
184     }
185   } 
186   delete [] index;
187   delete [] indexMin;
188   delete [] indexMax;
189   if(valden>0)val=valnum/valden;
190   AliInfo(Form(" the Average Efficiency = %f ",val)); 
191   return val;
192
193 //____________________________________________________________________
194 void AliCFEffGrid::Copy(TObject& eff) const
195 {
196   //
197   // copy function
198   //
199   Copy(eff);
200   AliCFEffGrid& target = (AliCFEffGrid &) eff;
201   
202   target.fSelNum=fSelNum; 
203   target.fSelDen=fSelDen; 
204   if(fContainer)
205     target.fContainer=fContainer;
206 }
207 //___________________________________________________________________
208 TH1F *AliCFEffGrid::Project(Int_t ivar) const
209 {
210   //
211   // Make a 1D projection along variable ivar 
212   //
213  
214   TH1F *proj1D=0;
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]];
219   }
220   
221   char pname[30];
222   sprintf(pname,"%s%s%i%i%s%i",GetName(),"_SelStep",fSelNum,fSelDen,"_proj1D_var", ivar);
223   char htitle[30];
224   sprintf(htitle,"%s%s%i%i%s%i",GetName(),"_SelStep",fSelNum,fSelDen,"_proj1D_var", ivar);
225   
226   if(!proj1D){
227     proj1D =new TH1F(pname,htitle, nbins, bins);
228   }  
229   
230   proj1D->Divide(fContainer->GetGrid(fSelNum)->Project(ivar),fContainer->GetGrid(fSelDen)->Project(ivar),1.,1.,"B");
231   
232   delete [] bins; 
233   return proj1D;
234
235 //___________________________________________________________________
236 TH2F *AliCFEffGrid::Project(Int_t ivar1,Int_t ivar2) const
237 {
238   //
239   // Make a 2D projection along variable ivar1,ivar2 
240   //
241  
242   TH2F *proj2D=0;
243
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]];
250   }
251   for(Int_t ibin2 =0;ibin2<nbins2+1;ibin2++){
252     bins2[ibin2] = fVarBinLimits[ibin2+fOffset[ivar2]];
253   }
254   
255   char pname[30];
256   sprintf(pname,"%s%s%i%i%s%i%i",GetName(),"_SelStep",fSelNum,fSelDen,"_proj2D_var", ivar1,ivar2);
257   char htitle[30];
258   sprintf(htitle,"%s%s%i%i%s%i%i",GetName(),"_SelStep",fSelNum,fSelDen,"_proj2D_var",ivar1,ivar2);
259   
260   if(!proj2D){
261     proj2D =new TH2F(pname,htitle, nbins1,bins1,nbins2,bins2);
262   }  
263   
264   proj2D->Divide(fContainer->GetGrid(fSelNum)->Project(ivar1,ivar2),fContainer->GetGrid(fSelDen)->Project(ivar1,ivar2),1.,1.,"B");
265   
266   delete [] bins1;
267   delete [] bins2; 
268   return proj2D;
269
270 //___________________________________________________________________
271 TH3F *AliCFEffGrid::Project(Int_t ivar1, Int_t ivar2, Int_t ivar3) const
272 {
273   //
274   // Make a 3D projection along variable ivar1,ivar2,ivar3 
275   //
276
277   TH3F *proj3D=0;
278
279   Int_t nbins1 =fNVarBins[ivar1];
280   Int_t nbins2 =fNVarBins[ivar2];
281   Int_t nbins3 =fNVarBins[ivar3];
282   
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];     
286   
287   for(Int_t ibin =0;ibin<nbins1+1;ibin++){
288     bins1[ibin] = fVarBinLimits[ibin+fOffset[ivar1]];
289   }
290   for(Int_t ibin =0;ibin<nbins2+1;ibin++){
291     bins2[ibin] = fVarBinLimits[ibin+fOffset[ivar2]];
292   }
293   for(Int_t ibin =0;ibin<nbins3+1;ibin++){
294     bins3[ibin] = fVarBinLimits[ibin+fOffset[ivar3]];
295   }
296   
297   char pname[30];
298   sprintf(pname,"%s%s%i%i%s%i%i%i",GetName(),"_SelStep",fSelNum,fSelDen,"_proj3D_var",ivar1,ivar2,ivar3);
299   char htitle[30];
300   sprintf(htitle,"%s%s%i%i%s%i%i%i",GetName(),"_SelStep",fSelNum,fSelDen,"_proj3D_var",ivar1,ivar2,ivar3);
301    
302   
303   if(!proj3D){
304     proj3D =new TH3F(pname,htitle, nbins1, bins1,nbins2,bins2,nbins3,bins3);
305   }  
306   
307   proj3D->Divide(fContainer->GetGrid(fSelNum)->Project(ivar1,ivar2,ivar3),fContainer->GetGrid(fSelDen)->Project(ivar1,ivar2,ivar3),1.,1.,"B");
308   
309   delete [] bins1;
310   delete [] bins2;
311   delete [] bins3; 
312   
313   return proj3D;
314