]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/AliCFEffGrid.cxx
Revert unwanted changes concerning PID in HLT
[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 "TMath.h"
29 #include "AliLog.h"
30 #include "AliCFEffGrid.h"
31 #include "TH1D.h"
32 #include "TH2D.h"
33 #include "TH3D.h"
34
35 //____________________________________________________________________
36 ClassImp(AliCFEffGrid)
37
38 //____________________________________________________________________
39 AliCFEffGrid::AliCFEffGrid() : 
40   AliCFGridSparse(),
41   fContainer(0x0),
42   fSelNum(-1),
43   fSelDen(-1)
44 {
45   //
46   // default constructor
47   //
48 }
49
50 //____________________________________________________________________
51 AliCFEffGrid::AliCFEffGrid(const Char_t* name, const Char_t* title, const Int_t nVarIn, const Int_t * nBinIn) :
52   AliCFGridSparse(name,title,nVarIn,nBinIn),
53   fContainer(0x0),
54   fSelNum(-1),
55   fSelDen(-1)
56 {
57   //
58   // ctor
59   //
60   SumW2();
61 }
62 //____________________________________________________________________
63 AliCFEffGrid::AliCFEffGrid(const Char_t* name, const Char_t* title, const AliCFContainer &c) :  
64   AliCFGridSparse(name,title,c.GetNVar(),c.GetNBins()),
65   fContainer(NULL),
66   fSelNum(-1),
67   fSelDen(-1)
68 {
69   //
70   // main constructor
71   //
72   SumW2();
73   //assign the container;
74   fContainer=&c;
75   for (Int_t iVar=0; iVar<GetNVar(); iVar++) {
76     SetBinLimits(iVar,GetBinLimits(iVar));
77   }
78   for (Int_t iVar=0; iVar<GetNVar(); iVar++) SetVarTitle(iVar,c.GetVarTitle(iVar));
79 }
80 //____________________________________________________________________
81 AliCFEffGrid::AliCFEffGrid(const AliCFEffGrid& eff) : 
82   AliCFGridSparse(eff),
83   fContainer(0x0),
84   fSelNum(-1),
85   fSelDen(-1)
86 {
87   //
88   // copy constructor
89   //
90   ((AliCFEffGrid &)eff).Copy(*this);
91 }
92
93 //____________________________________________________________________
94 AliCFEffGrid::~AliCFEffGrid()
95 {
96   //
97   // destructor
98   //
99 }
100
101 //____________________________________________________________________
102 AliCFEffGrid &AliCFEffGrid::operator=(const AliCFEffGrid &eff)
103 {
104   //
105   // assigment operator
106   //
107   if (this != &eff) eff.Copy(*this);
108   return *this;
109
110
111 //____________________________________________________________________
112 void AliCFEffGrid::CalculateEfficiency(Int_t istep1,Int_t istep2, Option_t *option)
113 {
114   //
115   // Calculate the efficiency matrix and its error between selection
116   // Steps istep1 and istep2
117   //
118   // 'option' is used as an argument for THnSparse::Divide
119   // default is "B" : binomial error calculation
120   //
121
122   fSelNum=istep1;
123   fSelDen=istep2;
124   AliCFGridSparse *num=GetNum();
125   AliCFGridSparse *den=GetDen();
126   num->SumW2();
127   den->SumW2();
128   this->SumW2();
129   this->Divide(num,den,1.,1.,option);
130   SetTitle(Form("Efficiency: %s / %s",fContainer->GetStepTitle(istep1),fContainer->GetStepTitle(istep2)));
131
132   AliInfo(Form("Efficiency calculated for steps %i and %i.",fSelNum,fSelDen));
133
134 //_____________________________________________________________________
135 Double_t AliCFEffGrid::GetAverage() const 
136 {
137   //
138   // Get the average efficiency 
139   //
140
141   Double_t val=0;
142   Double_t valnum=0;
143   Double_t valden=0;
144
145   THnSparse* num = ((AliCFGridSparse*)GetNum())->GetGrid() ;
146   THnSparse* den = ((AliCFGridSparse*)GetDen())->GetGrid() ;
147
148   for (Long_t iBin=0; iBin<num->GetNbins(); iBin++) valnum+=num->GetBinContent(iBin);
149   for (Long_t iBin=0; iBin<den->GetNbins(); iBin++) valden+=den->GetBinContent(iBin);
150   if (valden>0) val=valnum/valden;
151   AliInfo(Form(" The Average Efficiency = %f ",val)); 
152   return val;
153
154 //_____________________________________________________________________
155 // Double_t AliCFEffGrid::GetAverage(const Double_t *varMin, const Double_t* varMax ) const 
156 // {
157 //   //
158 //   // Get ave efficiency in a range
159 //   // (may not work properly, should be modified)
160
161 //   Double_t val=0;
162 //   Int_t *indexMin = new Int_t[GetNVar()];
163 //   Int_t *indexMax = new Int_t[GetNVar()];
164 //   Int_t *index    = new Int_t[GetNVar()];
165   
166 //   //Find out the min and max bins
167   
168 //   for(Int_t i=0;i<GetNVar();i++){
169 //     Double_t xmin=varMin[i]; // the min values  
170 //     Double_t xmax=varMax[i]; // the max values  
171 //     Int_t nbins = GetNBins(i)+1;
172 // //     Double_t *bins = new Double_t[nbins];
173 // //     for (Int_t ibin =0; ibin<nbins; ibin++) {
174 // //       bins[ibin] = fVarBinLimits[ibin+fOffset[i]];
175 // //     }
176 //     bins = GetBinLimits(i);
177 //     indexMin[i] = TMath::BinarySearch(nbins,bins,xmin);
178 //     indexMax[i] = TMath::BinarySearch(nbins,bins,xmax);
179 //     if(xmax>=bins[nbins-1]){
180 //       indexMax[i]=indexMax[i]-1;
181 //     }  
182 //     delete [] bins;
183 //   }
184   
185 //   Double_t valnum=0;
186 //   Double_t valden=0;
187 //   for (Int_t i=0; i<fNDim; i++) {
188 //     for (Int_t j=0; j<GetNVar(); j++) index[j]=GetBinIndex(j,i);
189 //     Bool_t isIn=kTRUE;
190 //     for (Int_t j=0;j<GetNVar();j++) {
191 //       if(!(index[j]>=indexMin[j] && index[j]<=indexMax[j]))isIn=kFALSE;   
192 //     }
193 //     if(isIn){
194 //       valnum+=GetNum()->GetElement(i);
195 //       valden+=GetDen()->GetElement(i);
196 //     }
197 //   } 
198 //   delete [] index;
199 //   delete [] indexMin;
200 //   delete [] indexMax;
201 //   if(valden>0)val=valnum/valden;
202 //   AliInfo(Form(" the Average Efficiency = %f ",val)); 
203 //   return val;
204 // } 
205 //____________________________________________________________________
206 // void AliCFEffGrid::Copy(TObject& eff) const
207 // {
208 //   //
209 //   // copy function
210 //   //
211 //   Copy(eff);
212 //   AliCFEffGrid& target = (AliCFEffGrid &) eff;
213   
214 //   target.fSelNum=fSelNum; 
215 //   target.fSelDen=fSelDen; 
216 //   if(fContainer)
217 //     target.fContainer=fContainer;
218 // }
219 //___________________________________________________________________
220 TH1D *AliCFEffGrid::Project(Int_t ivar) const
221 {
222   //
223   // Make a 1D projection along variable ivar 
224   //
225   
226   if (fSelNum<0 || fSelDen<0) {
227     AliError("You must call CalculateEfficiency() first !");
228     return 0x0;
229   }
230   const Int_t nDim = 1 ;
231   Int_t dim[nDim] = {ivar} ;
232   THnSparse* hNum = ((AliCFGridSparse*)GetNum())->GetGrid()->Projection(nDim,dim);
233   THnSparse* hDen = ((AliCFGridSparse*)GetDen())->GetGrid()->Projection(nDim,dim);
234   THnSparse* ratio = (THnSparse*)hNum->Clone();
235   ratio->Divide(hNum,hDen,1.,1.,"B");
236   delete hNum; delete hDen;
237   TH1D* h = ratio->Projection(0);
238   h->SetXTitle(GetVarTitle(ivar));
239   h->SetTitle(Form("%s projected on %s",GetTitle(),GetVarTitle(ivar)));
240   return h ;
241
242 //___________________________________________________________________
243 TH2D *AliCFEffGrid::Project(Int_t ivar1,Int_t ivar2) const
244 {
245   //
246   // Make a 2D projection along variable ivar1,ivar2 
247   //
248   
249   if (fSelNum<0 || fSelDen<0) {
250     AliError("You must call CalculateEfficiency() first !");
251     return 0x0;
252   }
253   const Int_t nDim = 2 ;
254   Int_t dim[nDim] = {ivar1,ivar2} ;
255   THnSparse* hNum = ((AliCFGridSparse*)GetNum())->GetGrid()->Projection(nDim,dim);
256   THnSparse* hDen = ((AliCFGridSparse*)GetDen())->GetGrid()->Projection(nDim,dim);
257   THnSparse* ratio = (THnSparse*)hNum->Clone();
258   ratio->Divide(hNum,hDen,1.,1.,"B");
259   delete hNum; delete hDen;
260   TH2D* h = ratio->Projection(0,1);
261   h->SetXTitle(GetVarTitle(ivar1));
262   h->SetYTitle(GetVarTitle(ivar2));
263   h->SetTitle(Form("%s projected on %s-%s",GetTitle(),GetVarTitle(ivar1),GetVarTitle(ivar2)));
264   return h;
265
266 //___________________________________________________________________
267 TH3D *AliCFEffGrid::Project(Int_t ivar1, Int_t ivar2, Int_t ivar3) const
268 {
269   //
270   // Make a 3D projection along variable ivar1,ivar2,ivar3 
271   //
272
273   if (fSelNum<0 || fSelDen<0) {
274     AliError("You must call CalculateEfficiency() first !");
275     return 0x0;
276   }
277   const Int_t nDim = 3 ;
278   Int_t dim[nDim] = {ivar1,ivar2,ivar3} ;
279   THnSparse* hNum = ((AliCFGridSparse*)GetNum())->GetGrid()->Projection(nDim,dim);
280   THnSparse* hDen = ((AliCFGridSparse*)GetDen())->GetGrid()->Projection(nDim,dim);
281   THnSparse* ratio = (THnSparse*)hNum->Clone();
282   ratio->Divide(hNum,hDen,1.,1.,"B");
283   delete hNum; delete hDen;
284   TH3D* h = ratio->Projection(0,1,2);
285   h->SetXTitle(GetVarTitle(ivar1));
286   h->SetYTitle(GetVarTitle(ivar2));
287   h->SetZTitle(GetVarTitle(ivar3));
288   h->SetTitle(Form("%s projected on %s-%s-%s",GetTitle(),GetVarTitle(ivar1),GetVarTitle(ivar2),GetVarTitle(ivar3)));
289   return h;
290
291 //___________________________________________________________________
292 AliCFEffGrid* AliCFEffGrid::MakeSlice(Int_t nVars, const Int_t* vars, const Double_t* varMin, const Double_t* varMax, Int_t numStep, Int_t denStep) const {
293   //
294   // Makes a slice along the "nVars" variables defined in the array "vars[nVars]" for all the container steps.
295   // The ranges of ALL the container variables must be defined in the array varMin[fNVar] and varMax[fNVar].
296   // This function returns the effiency relative to this new 'sliced' container, between steps defined in numStep and denStep
297   //
298   
299   AliCFContainer* cont = fContainer->MakeSlice(nVars,vars,varMin,varMax);
300   AliCFEffGrid  * eff  = new AliCFEffGrid(Form("%s_sliced",GetName()), Form("%s_sliced",GetTitle()), *cont);
301   eff->CalculateEfficiency(numStep,denStep);
302   return eff;
303 }