]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/AliCFEffGrid.cxx
Fixed Coding Conventions for class AliACORDEv1 (Mario RC)
[u/mrichter/AliRoot.git] / CORRFW / AliCFEffGrid.cxx
1 /* $Id$ */
2 /**************************************************************************
3  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  *                                                                        *
5  * Author: The ALICE Off-line Project.                                    *
6  * Contributors are mentioned in the code where appropriate.              *
7  *                                                                        *
8  * Permission to use, copy, modify and distribute this software and its   *
9  * documentation strictly for non-commercial purposes is hereby granted   *
10  * without fee, provided that the above copyright notice appears in all   *
11  * copies and that both the copyright notice and this permission notice   *
12  * appear in the supporting documentation. The authors make no claims     *
13  * about the suitability of this software for any purpose. It is          *
14  * provided "as is" without express or implied warranty.                  *
15  **************************************************************************/
16 //--------------------------------------------------------------------//
17 //                                                                    //
18 // AliCFEffGrid Class                                              //
19 // Class to handle efficiency grids                                   // 
20 //                                                                    //
21 // -- Author : S.Arcelli                                              //
22 //                                                                    //
23 //                                                                    //
24 //                                                                    //
25 //--------------------------------------------------------------------//
26 //
27 //
28 #include <TROOT.h>
29 #include <TMath.h>
30 #include <TFile.h>
31 #include <AliLog.h>
32 #include "AliCFEffGrid.h"
33 #include <TH1D.h>
34 #include <TH2D.h>
35 #include <TH3D.h>
36
37 //____________________________________________________________________
38 ClassImp(AliCFEffGrid)
39
40 //____________________________________________________________________
41 AliCFEffGrid::AliCFEffGrid() : 
42   AliCFGrid(),
43   fContainer(0x0),
44   fSelNum(-1),
45   fSelDen(-1)
46 {
47   //
48   // default constructor
49   //
50 }
51
52 //____________________________________________________________________
53 AliCFEffGrid::AliCFEffGrid(const Char_t* name, const Char_t* title, const Int_t nVarIn, const Int_t * nBinIn, const Double_t *binLimitsIn) :  
54   AliCFGrid(name,title,nVarIn,nBinIn,binLimitsIn),
55   fContainer(0x0),
56   fSelNum(-1),
57   fSelDen(-1)
58 {
59   //
60   // ctor
61   //
62   SumW2();
63 }
64 //____________________________________________________________________
65 AliCFEffGrid::AliCFEffGrid(const Char_t* name, const Char_t* title, const AliCFContainer &c) :  
66   AliCFGrid(name,title,c.GetNVar(),c.GetNBins(),c.GetBinLimits()),
67   fContainer(0x0),
68   fSelNum(-1),
69   fSelDen(-1)
70 {
71   //
72   // main constructor
73   //
74   SumW2();
75   //assign the container;
76   fContainer=&c;
77 }
78 //____________________________________________________________________
79 AliCFEffGrid::AliCFEffGrid(const AliCFEffGrid& eff) :   AliCFGrid(),
80   fContainer(0x0),
81   fSelNum(-1),
82   fSelDen(-1)
83 {
84   //
85   // copy constructor
86   //
87   ((AliCFEffGrid &)eff).Copy(*this);
88 }
89
90 //____________________________________________________________________
91 AliCFEffGrid::~AliCFEffGrid()
92 {
93   //
94   // destructor
95   //
96 }
97
98 //____________________________________________________________________
99 AliCFEffGrid &AliCFEffGrid::operator=(const AliCFEffGrid &eff)
100 {
101   //
102   // assigment operator
103   //
104   if (this != &eff)
105     ((AliCFEffGrid &) eff).Copy(*this);
106   return *this;
107
108 //____________________________________________________________________
109
110 void AliCFEffGrid::CalculateEfficiency(Int_t istep1,Int_t istep2)
111 {
112   //
113   // Calculate the efficiency matrix and its error between selection
114   // Steps istep1 and istep2
115   //
116
117   fSelNum=istep1;
118   fSelDen=istep2;
119   AliCFVGrid *num=fContainer->GetGrid(fSelNum);
120   AliCFVGrid *den=fContainer->GetGrid(fSelDen);
121   num->SumW2();
122   den->SumW2();
123   this->SumW2();
124   this->Divide(num,den,1.,1.,"B");
125
126   Int_t nEmptyBinsNum=0;
127   Int_t nEmptyBinsNumAndDen=0;
128   for(Int_t iel=0;iel<fNDim;iel++){
129     if(den->GetElement(iel)>0){
130       if(num->GetElement(iel)==0)nEmptyBinsNum++; //num==0,den!=0
131     }
132     else{
133       nEmptyBinsNumAndDen++;
134     }
135   }    
136   // Some monitoring printout:
137   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));
138   AliInfo(Form("The correction map contains %i empty bins ",nEmptyBinsNum+nEmptyBinsNumAndDen));
139
140 //_____________________________________________________________________
141 Double_t AliCFEffGrid::GetAverage() const 
142 {
143   //
144   // Get the average efficiency 
145   //
146
147   Double_t val=0;
148   Double_t valnum=0;
149   Double_t valden=0;
150   for(Int_t i=0;i<fNDim;i++){
151     valnum+=fContainer->GetGrid(fSelNum)->GetElement(i);
152     valden+=fContainer->GetGrid(fSelDen)->GetElement(i);
153   }
154   if(valden>0)val=valnum/valden;
155   AliInfo(Form(" The Average Efficiency = %f ",val)); 
156
157   return val;
158
159 //_____________________________________________________________________
160 Double_t AliCFEffGrid::GetAverage(Double_t *varMin, Double_t* varMax ) const 
161 {
162   //
163   // Get ave efficiency in a range
164   //
165
166
167   Double_t val=0;
168   Int_t *indexMin = new Int_t[fNVar];
169   Int_t *indexMax = new Int_t[fNVar];
170   Int_t *index    = new Int_t[fNVar];
171   
172   //Find out the min and max bins
173   
174   for(Int_t i=0;i<fNVar;i++){
175     Double_t xmin=varMin[i]; // the min values  
176     Double_t xmax=varMax[i]; // the max values  
177     Int_t nbins=fNVarBins[i]+1;
178     Float_t *bins=new Float_t[nbins];
179     for(Int_t ibin =0;ibin<nbins;ibin++){
180       bins[ibin] = fVarBinLimits[ibin+fOffset[i]];
181     }
182     indexMin[i] = TMath::BinarySearch(nbins,bins,xmin);
183     indexMax[i] = TMath::BinarySearch(nbins,bins,xmax);
184     if(xmax>=bins[nbins-1]){
185       indexMax[i]=indexMax[i]-1;
186     }  
187     delete [] bins;
188   }
189   
190   Double_t valnum=0;
191   Double_t valden=0;
192   for(Int_t i=0;i<fNDim;i++){
193     for (Int_t j=0;j<fNVar;j++)index[j]=GetBinIndex(j,i);
194     Bool_t isIn=kTRUE;
195     for (Int_t j=0;j<fNVar;j++){
196       if(!(index[j]>=indexMin[j] && index[j]<=indexMax[j]))isIn=kFALSE;   
197     }
198     if(isIn){
199       valnum+=fContainer->GetGrid(fSelNum)->GetElement(i);
200       valden+=fContainer->GetGrid(fSelDen)->GetElement(i);
201     }
202   } 
203   delete [] index;
204   delete [] indexMin;
205   delete [] indexMax;
206   if(valden>0)val=valnum/valden;
207   AliInfo(Form(" the Average Efficiency = %f ",val)); 
208   return val;
209
210 //____________________________________________________________________
211 void AliCFEffGrid::Copy(TObject& eff) const
212 {
213   //
214   // copy function
215   //
216   Copy(eff);
217   AliCFEffGrid& target = (AliCFEffGrid &) eff;
218   
219   target.fSelNum=fSelNum; 
220   target.fSelDen=fSelDen; 
221   if(fContainer)
222     target.fContainer=fContainer;
223 }
224 //___________________________________________________________________
225 TH1D *AliCFEffGrid::Project(Int_t ivar) const
226 {
227   //
228   // Make a 1D projection along variable ivar 
229   //
230  
231   TH1D *proj1D=0;
232   Int_t nbins =fNVarBins[ivar];
233   Float_t *bins = new Float_t[nbins+1];    
234   for(Int_t ibin =0;ibin<nbins+1;ibin++){
235     bins[ibin] = fVarBinLimits[ibin+fOffset[ivar]];
236   }
237   
238   char pname[30];
239   sprintf(pname,"%s%s%i%i%s%i",GetName(),"_SelStep",fSelNum,fSelDen,"_proj1D_var", ivar);
240   char htitle[30];
241   sprintf(htitle,"%s%s%i%i%s%i",GetName(),"_SelStep",fSelNum,fSelDen,"_proj1D_var", ivar);
242   
243   if(!proj1D){
244     proj1D =new TH1D(pname,htitle, nbins, bins);
245   }  
246   
247   proj1D->Sumw2();
248   proj1D->Divide(fContainer->GetGrid(fSelNum)->Project(ivar),fContainer->GetGrid(fSelDen)->Project(ivar),1.,1.,"B");
249   
250   delete [] bins; 
251   return proj1D;
252
253 //___________________________________________________________________
254 TH2D *AliCFEffGrid::Project(Int_t ivar1,Int_t ivar2) const
255 {
256   //
257   // Make a 2D projection along variable ivar1,ivar2 
258   //
259  
260   TH2D *proj2D=0;
261
262   Int_t nbins1 =fNVarBins[ivar1];
263   Float_t *bins1 = new Float_t[nbins1+1];    
264   Int_t nbins2 =fNVarBins[ivar2];
265   Float_t *bins2 = new Float_t[nbins2+1];    
266   for(Int_t ibin1 =0;ibin1<nbins1+1;ibin1++){
267     bins1[ibin1] = fVarBinLimits[ibin1+fOffset[ivar1]];
268   }
269   for(Int_t ibin2 =0;ibin2<nbins2+1;ibin2++){
270     bins2[ibin2] = fVarBinLimits[ibin2+fOffset[ivar2]];
271   }
272   
273   char pname[30];
274   sprintf(pname,"%s%s%i%i%s%i%i",GetName(),"_SelStep",fSelNum,fSelDen,"_proj2D_var", ivar1,ivar2);
275   char htitle[30];
276   sprintf(htitle,"%s%s%i%i%s%i%i",GetName(),"_SelStep",fSelNum,fSelDen,"_proj2D_var",ivar1,ivar2);
277   
278   if(!proj2D){
279     proj2D =new TH2D(pname,htitle, nbins1,bins1,nbins2,bins2);
280   }  
281   
282   proj2D->Sumw2();
283   proj2D->Divide(fContainer->GetGrid(fSelNum)->Project(ivar1,ivar2),fContainer->GetGrid(fSelDen)->Project(ivar1,ivar2),1.,1.,"B");
284   
285   delete [] bins1;
286   delete [] bins2; 
287   return proj2D;
288
289 //___________________________________________________________________
290 TH3D *AliCFEffGrid::Project(Int_t ivar1, Int_t ivar2, Int_t ivar3) const
291 {
292   //
293   // Make a 3D projection along variable ivar1,ivar2,ivar3 
294   //
295
296   TH3D *proj3D=0;
297
298   Int_t nbins1 =fNVarBins[ivar1];
299   Int_t nbins2 =fNVarBins[ivar2];
300   Int_t nbins3 =fNVarBins[ivar3];
301   
302   Float_t *bins1 = new Float_t[nbins1+1];         
303   Float_t *bins2 = new Float_t[nbins2+1];     
304   Float_t *bins3 = new Float_t[nbins3+1];     
305   
306   for(Int_t ibin =0;ibin<nbins1+1;ibin++){
307     bins1[ibin] = fVarBinLimits[ibin+fOffset[ivar1]];
308   }
309   for(Int_t ibin =0;ibin<nbins2+1;ibin++){
310     bins2[ibin] = fVarBinLimits[ibin+fOffset[ivar2]];
311   }
312   for(Int_t ibin =0;ibin<nbins3+1;ibin++){
313     bins3[ibin] = fVarBinLimits[ibin+fOffset[ivar3]];
314   }
315   
316   char pname[30];
317   sprintf(pname,"%s%s%i%i%s%i%i%i",GetName(),"_SelStep",fSelNum,fSelDen,"_proj3D_var",ivar1,ivar2,ivar3);
318   char htitle[30];
319   sprintf(htitle,"%s%s%i%i%s%i%i%i",GetName(),"_SelStep",fSelNum,fSelDen,"_proj3D_var",ivar1,ivar2,ivar3);
320    
321   
322   if(!proj3D){
323     proj3D =new TH3D(pname,htitle, nbins1, bins1,nbins2,bins2,nbins3,bins3);
324   }  
325   
326   proj3D->Sumw2();
327   proj3D->Divide(fContainer->GetGrid(fSelNum)->Project(ivar1,ivar2,ivar3),fContainer->GetGrid(fSelDen)->Project(ivar1,ivar2,ivar3),1.,1.,"B");
328   
329   delete [] bins1;
330   delete [] bins2;
331   delete [] bins3; 
332   
333   return proj3D;
334