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