]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FASTSIM/AliFastMuonTriggerEff.cxx
Adding MUON HLT code to the repository.
[u/mrichter/AliRoot.git] / FASTSIM / AliFastMuonTriggerEff.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 #include "AliFastMuonTriggerEff.h"
19 #include <TFile.h>
20 #include <TTree.h>
21 #include <TROOT.h>
22 #include <stdlib.h>
23 #include <Rtypes.h>
24
25 // Debugging flag
26 //#define MYTRIGDEBUG
27
28 ClassImp(AliFastMuonTriggerEff)
29
30 AliFastMuonTriggerEff::AliFastMuonTriggerEff():
31     AliFastResponse("Efficiency", "Muon Trigger Efficiency")
32 {
33 //
34 // Default constructor
35 //
36     SetCut(kLow);
37     fZones=0;
38 }
39
40
41 void AliFastMuonTriggerEff::Init()
42 {
43 //
44 //  Initialization
45 //
46     printf("Initializing %s / %s", fName.Data(), fTitle.Data());
47 //
48
49     fPhiMin   = -90. ;
50     fPhiMax   =  90. ;
51     fDphi     =   9. ;
52     fThetaMin =   2. ;
53     fThetaMax =   9. ;
54     fDtheta   =   0.7;
55     fDpt      =   0.1;
56
57     InitTree();
58 }
59
60 void AliFastMuonTriggerEff::InitTree()
61 {
62 //
63 //  Initialize tables from tree
64 //
65     TTree          *chain;   // pointer to the analyzed TTree or TChain
66     Int_t           nSim, nZona;
67     Double_t        pmin, pmax, tmin, tmax, phimin, phimax, bkg;
68     Double_t        len50, hen50, leff, heff;
69     Double_t        vLPt[50];
70     Double_t        vHPt[50];
71     Char_t file[100]="$(ALICE_ROOT)/FASTSIM/data/vettorpara.root";
72 //
73 //  Avoid memory leak in case of reinitialization
74     if(fZones!=0) {
75         printf("\nWarning: reinitialization of an object of class: AliFastMuonTriggerEff\n");
76         for (Int_t i=0; i<fZones; i++) {
77             if(fEffLow [i])delete[] fEffLow[i];
78             if(fEffHigh[i])delete[] fEffHigh[i];
79         }
80         if(fEffLow) {
81             delete[] fEffLow;
82             fEffLow=0;
83         }
84         if(fEffHigh) {
85             delete[] fEffHigh;
86             fEffHigh=0;
87         }
88     }
89     printf("AliFastMuonTriggerEff: Initialization\n");
90     TFile *f = new TFile(file);
91     if(f->IsZombie()) {
92         printf("Cannot open file: %s\n",file);
93         return;
94     }
95     f->ls();
96     
97     TTree* tree = (TTree*)gDirectory->Get("fitPar");
98
99 //
100 //
101     chain = tree;
102     chain->SetMakeClass(1);
103     
104     chain->SetBranchAddress("nSim",&nSim);
105     chain->SetBranchAddress("nZona",&nZona);
106     chain->SetBranchAddress("ptmin",&pmin);
107     chain->SetBranchAddress("ptmax",&pmax);
108     chain->SetBranchAddress("Thetamin",&tmin);
109     chain->SetBranchAddress("Thetamax",&tmax);
110     chain->SetBranchAddress("Phimin",&phimin);
111     chain->SetBranchAddress("Phimax",&phimax);
112     chain->SetBranchAddress("Bkg",&bkg);
113     chain->SetBranchAddress("EffLPt",vLPt);
114     chain->SetBranchAddress("EffHPt",vHPt);
115     chain->SetBranchAddress("Pt0.5Low",&len50);
116     chain->SetBranchAddress("Pt0.5High",&hen50);
117     chain->SetBranchAddress("EffLowMax",&leff);
118     chain->SetBranchAddress("EffHighMax",&heff);
119 // 
120 //
121     Int_t nentries = Int_t(chain->GetEntries());
122     Int_t nbytes = 0, nb = 0;
123
124 // Count the number of zones of the parametrization
125     Int_t nzone0=0, nzone1=0;
126     for (Int_t jentry=0; jentry<nentries; jentry++) {
127         nb = chain->GetEntry(jentry);
128         if(bkg==0.) {
129              if(nSim==0)nzone0++;
130              if(nSim==1)nzone1++;
131         }
132     }
133     
134     printf("Trigger parametrisation for %d zones for Pt: 0. 5. GeV/c\n",nzone0);
135     printf("and %d zones extended to 10. GeV/c\n",nzone1);
136     fZones=nzone0+nzone1;
137 //    printf("Ciao\n");
138     if(fZones<=0){
139         printf("Number of zones must be greater than 0\n");
140         exit(6);
141     }
142     
143     fEffLow =new Float_t*[fZones];
144     fEffHigh=new Float_t*[fZones];
145     for (Int_t i=0; i<fZones; i++) {
146         fEffLow [i]=new Float_t[50];
147         fEffHigh[i]=new Float_t[50];
148     }
149     
150 //  Initialize look-up table to standard values
151     Int_t isim, itheta, iphi;
152     for (isim=0; isim<2; isim++) {
153         for (itheta=0; itheta<10; itheta++) {
154             for (iphi=0; iphi<20; iphi++) {
155                 fLook[isim][itheta][iphi]=0;
156             }
157         }
158     }
159
160 //  Loading Trigger efficiencies
161     Int_t myzone=0;
162 #ifdef MYTRIGDEBUG
163             printf("Filling nSim nZona pmin pmax tmin tmax phimin phimax: ....\n");
164 #endif
165     for (Int_t jentry=0; jentry<nentries; jentry++) {
166 //      Int_t ientry = LoadTree(jentry); 
167         nb = chain->GetEntry(jentry);   
168         nbytes += nb;
169 #ifdef MYTRIGDEBUG
170        printf("Getting entry %d... ",jentry);
171 #endif
172 // For the time being it works with background 0
173         if ((nSim == 0 || nSim == 1)&&bkg==0.) {
174 #ifdef MYTRIGDEBUG
175             printf("Filling %d %d %6.2f %6.2f %6.2f %6.2f %6.2f %6.2f: ",
176             nSim,nZona,pmin,pmax,tmin,tmax,phimin,phimax);
177 #endif
178             for (Int_t k = 0; k < 50; k++) {
179                 fEffLow [myzone][k] = vLPt[k];
180                 fEffHigh[myzone][k] = vHPt[k];
181 #ifdef MYTRIGDEBUG
182                 if(k<15)printf(" %5.3f",vLPt[k]);
183 #endif
184             }
185 #ifdef MYTRIGDEBUG
186             printf("\n");
187 #endif
188             myzone++;
189             iphi=Int_t(nZona/10)-10;
190             itheta=nZona%10;
191             if(iphi<0||iphi>19||itheta<0||itheta>9) {
192                 printf("The format of data file is not consistent\nwith this version of the code\n");
193                 printf("This should never happen: iphi %d, itheta: %d\n",iphi,itheta);
194                 exit(7);
195             }
196             fLook[nSim][itheta][iphi]=myzone;
197         } else {
198             printf("Skipping entry with nSim=%d, bkg=%f\n",nSim,bkg);
199         }
200     }
201 #ifdef MYTRIGDEBUG
202     printf("This is the content of the LUT after first step of initialization\n");
203     for(isim=0; isim<2; isim++) {
204         printf("isim=%d\n",isim);
205         for(iphi=0; iphi<20; iphi++) {
206             printf("iphi: %2d:",iphi);
207             for(itheta=0; itheta<10; itheta++) {
208                 printf(" %4d",fLook[isim][itheta][iphi]);
209             }
210             printf("\n");
211         }
212     }
213 #endif
214 // Filling look=up table for the zones where the extended simulation does
215 // not exists
216     for(iphi=0; iphi<20; iphi++) {
217         for(itheta=0; itheta<10; itheta++) {
218             if(fLook[0][itheta][iphi]==0) {
219                 printf("Missing entry isim=%d itheta=%d iphi=%d\n",isim,itheta,iphi);
220                 exit(8);
221             }
222             if(fLook[0][itheta][iphi]<0||fLook[0][itheta][iphi]>fZones) {
223                 printf("Problem with entry isim=%d itheta=%d iphi=%d\n",isim,itheta,iphi);
224                 exit(9);
225             }
226             if(fLook[1][itheta][iphi]==0) {
227                  fLook[1][itheta][iphi]=-fLook[0][itheta][iphi];
228             }
229         }
230     }
231 #ifdef MYTRIGDEBUG
232     for(isim=0; isim<2; isim++) {
233         printf("isim=%d\n",isim);
234         for(iphi=0; iphi<20; iphi++) {
235             printf("iphi: %2d:",iphi);
236             for(itheta=0; itheta<10; itheta++) {
237                 printf(" %4d",fLook[isim][itheta][iphi]);
238             }
239             printf("\n");
240         }
241     }
242         for(iphi=0; iphi<20; iphi++) {
243             for(itheta=0; itheta<10; itheta++) {
244     for(isim=0; isim<2; isim++) {
245                 printf("%d %2d %2d %4d:",isim,iphi,itheta,fLook[isim][itheta][iphi]);
246                 if(fLook[isim][itheta][iphi]>0) {
247                     myzone=fLook[isim][itheta][iphi]-1;
248                     for(Int_t ipt=0; ipt<20; ipt++) {
249                         //printf(" %5.3f",fEffLow[myzone][ipt]);
250                         printf(" %5.3f",fEffHigh[myzone][ipt]);
251                     }
252                     printf(" ...");
253                     for(Int_t ipt=40; ipt<50; ipt++) {
254                         //printf(" %5.3f",fEffLow[myzone][ipt]);
255                         printf(" %5.3f",fEffHigh[myzone][ipt]);
256                     }               
257                 }
258                 printf("\n");
259             }
260         }
261     }    
262 #endif
263     f->Close();
264 }
265
266 void AliFastMuonTriggerEff::Evaluate(Float_t charge, Float_t pt,
267                 Float_t theta, Float_t phi,Float_t& effLow, Float_t& effHigh,
268                 Float_t& /*eff*/)
269 {
270 //
271 //  Trigger efficiency for pt, theta, phi (low and high cut)
272 //
273     effLow=0.;
274     effHigh=0.;
275     if(fZones==0) {
276         printf("Call to uninitialized object of class: AliFastMuonTriggerEff\n");
277         return;
278     }
279     if(pt<0) {
280         printf("Warning: pt: %f < 0. GeV/c\n",pt);
281         return; 
282     }
283     Int_t iPt   = Int_t(pt/fDpt)%50;
284     Int_t iSim  = Int_t(pt/fDpt)/50;
285     Int_t iPhi  = Int_t((phi-fPhiMin)/fDphi);
286     if(phi<fPhiMin)iPhi=iPhi-1;
287     Int_t iTheta = Int_t((theta-fThetaMin)/fDtheta);
288 #ifdef MYTRIGDEBUG
289     printf("iSim iPt iTheta iPhi: %d %d %d %d\n",iSim,iPt,iTheta,iPhi);
290 #endif
291     iPhi=iPhi-40*(iPhi/40);
292 #ifdef MYTRIGDEBUG
293     printf("1: iPhi converted to: %d for angle equivalence\n",iPhi);
294 #endif
295     if(iPhi<0)iPhi=-iPhi-1;
296     if(iPhi>19)iPhi=39-iPhi;
297 #ifdef MYTRIGDEBUG
298     printf("2: iPhi converted to: %d for the symmetry of the spectrometer\n",iPhi);
299 #endif
300     if(charge==1.){
301     } else if(charge==-1.) {
302     iPhi=19-iPhi;
303 #ifdef MYTRIGDEBUG
304     printf("3: iPhi converted to: %d for charge symmetry\n",iPhi);
305 #endif
306     } else {
307         printf("Warning: not understand charge: %f\n",charge);
308         return;
309     }
310     if(iTheta<0||iTheta>9) {
311         printf("Warning: theta: %f outside acceptance\n",theta);
312         return;
313     }
314     if(iPt<0) {
315         printf("Warning: what do you mean with pt: %f <0?\n",pt);
316         return;
317     }
318     if(iSim>=fSim) {
319         iSim=fSim-1;
320         iPt=49;
321 #ifdef MYTRIGDEBUG
322     printf("4: iSim iPt converted to: %d %d (last zone)\n",iSim,iPt);
323 #endif
324     }
325     Int_t iLook=fLook[iSim][iTheta][iPhi];
326     if(iLook<0) {
327 #ifdef MYTRIGDEBUG
328     printf("5: iLook iPt: %d %d converted to: ",iLook,iPt);
329 #endif
330         iLook=-iLook-1;
331         iPt=49;
332 #ifdef MYTRIGDEBUG
333     printf("%d %d from look up table contents\n",iLook,iPt);
334 #endif
335     } else {
336         iLook=iLook-1;
337     }
338     effLow=fEffLow[iLook][iPt];
339     effHigh=fEffHigh[iLook][iPt];
340 #ifdef MYTRIGDEBUG
341     printf("6: charge, iSim, iTheta, iPhi, iPt: %f %d %d %d %d effLow: %f, effHigh: %f\n",
342                charge,iSim,iTheta,iPhi,iPt,effLow,effHigh);
343 #endif
344     
345     //fEffLow [iPhi][iTheta][iPt];
346     //fEffHigh[iPhi][iTheta][iPt];    
347 }
348
349 Float_t AliFastMuonTriggerEff::Evaluate(Float_t charge, Float_t pt,
350                    Float_t theta, Float_t phi)
351 {
352     //
353     // Trigger efficiency for pt, theta, phi depending of fCut
354     // 
355     if(fZones==0) {
356         printf("Call to uninitialized object of class: AliFastMuonTriggerEff\n");
357         return 0.;
358     }
359     Float_t eff;
360     Float_t effLow, effHigh;
361     
362     Evaluate(charge,pt,theta,phi,effLow,effHigh);
363     if (fCut == kLow) 
364         eff  = effLow;
365     else
366         eff  = effHigh;
367
368     return eff;
369 }