]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONTriggerEfficiencyPt.C
Protection against not accessible input OCDB (Laurent)
[u/mrichter/AliRoot.git] / MUON / MUONTriggerEfficiencyPt.C
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
19 // Macro to produce trigger single muon efficiency versus pt plots for the 
20 // 2 pt cuts. 
21 // see full description in the README file
22 // Author: Fabien Guerin (LPC)
23
24 // ROOT includes
25 #include "TBranch.h"
26 #include "TClonesArray.h"
27 #include "TFile.h"
28 #include "TF1.h"
29 #include "TH1.h"
30 #include "TH2.h"
31 #include "TParticle.h"
32 #include "TTree.h"
33 #include "TMath.h"
34 #include "TStyle.h"
35 #include "TCanvas.h"
36 #include "TLegend.h"
37 #include "Riostream.h"
38
39 // STEER includes
40 #include "AliRun.h"
41 #include "AliRunLoader.h"
42 #include "AliHeader.h"
43 #include "AliLoader.h"
44 #include "AliStack.h"
45
46 // MUON includes
47 #include "AliMUONDataInterface.h"
48 #include "AliMUONMCDataInterface.h"
49 #include "AliMUONVHitStore.h"
50 #include "AliMUONVTriggerStore.h"
51 #include "AliMUONHit.h"
52 #include "AliMUONDigit.h"
53 #include "AliMUONGlobalTrigger.h"
54 #include "AliMUONTrack.h"
55
56 Double_t fitArc(Double_t *x,Double_t *par) 
57 {
58   Double_t h1=par[1]*(x[0]+par[2]);
59   return par[0]*TMath::ATan(TMath::Power(h1,6))+par[3];
60 }
61
62 Double_t fitArch(Double_t *x,Double_t *par) 
63 {
64   Double_t h0=2*par[0]/(TMath::Pi());
65   Double_t h1=par[1]*x[0]+par[2]*x[0]*x[0];
66   return h0*TMath::TanH(h1)+par[3];
67 }
68
69 void MUONTriggerEfficiencyPt(const char* filenameSim="galice_sim.root", 
70                              const char* filenameRec="galice.root",  
71                              Bool_t readFromRP = 0)
72 {
73
74 // define style
75     TStyle *st1 = new TStyle("st1","My STYLE");
76     
77     st1->SetOptStat(0);
78     st1->SetOptFit(111);  
79     st1->SetFrameFillColor(10);
80     st1->SetCanvasColor(10);
81     st1->SetFillColor(10);
82     st1->SetTitleBorderSize(0);
83     st1->SetTitleOffset(1.2,"XY");
84     st1->SetTitleSize(0.05,"XY"); 
85     st1->SetLabelSize(0.045,"XY");  
86     st1->SetLabelFont(22,"XY");  
87     st1->SetTitleFont(22,"XY");     
88     st1->SetOptTitle(0);   
89     st1->SetStatColor(10);
90     st1->SetStatFont(62);     
91     st1->SetFitFormat("4.4g");
92     st1->SetPadLeftMargin(0.15);
93     st1->SetPadRightMargin(0.06); 
94     st1->SetPadTopMargin(0.06);
95     st1->SetPadBottomMargin(0.15); 
96     st1->cd();
97     
98 //    gROOT->ForceStyle();
99     //TGaxis::SetMaxDigits(3);
100     
101 // beginning of macro    
102     TH1F *ptlpt    = new TH1F("ptlpt","",50,0.15,5.); 
103     TH1F *pthpt    = new TH1F("pthpt","",50,0.15,5.);
104     TH1F *ptcoinc  = new TH1F("ptcoinc","",50,0.15,5.);  
105         
106     TParticle *particle;
107     
108     Int_t nevents;
109     Double_t coincmuon=0;
110     Double_t lptmuon=0;
111     Double_t hptmuon=0;
112     Double_t ptmu=0;
113
114     AliMUONMCDataInterface diSim(filenameSim);
115     AliMUONDataInterface diRec(filenameRec);
116     
117     if (!diSim.IsValid())
118     {
119         cerr << "Cannot access " << filenameSim << endl;
120         return;
121     }
122     
123     if (!diRec.IsValid())
124     {
125         cerr << "Cannot access " << filenameRec << endl;
126         return;
127     }
128     
129     FILE* fdat = fopen("MUONTriggerEfficiencyPt.out","w");          
130     if (!fdat)
131     {
132         cerr << "Cannot create output file" << endl;
133         return;
134     }
135
136     nevents = diSim.NumberOfEvents();
137
138     AliRunLoader* runLoader = AliRunLoader::Open(filenameSim);    
139
140     for (Int_t ievent=0; ievent<nevents; ievent++) {  // Event loop
141         
142         if (ievent%500==0) printf("ievent = %d \n",ievent);
143 // kine
144         
145         runLoader->GetRunLoader()->GetEvent(ievent);
146         runLoader->GetRunLoader()->LoadKinematics();
147         AliStack* stack = runLoader->GetRunLoader()->Stack();
148         
149         Int_t nparticles = (Int_t) stack->GetNtrack();        
150
151         for (Int_t iparticle=0; iparticle<nparticles; iparticle++) {
152             particle = stack->Particle(iparticle);   
153             Float_t pt = particle->Pt();            
154             Int_t pdgcode = particle->GetPdgCode();                        
155             if (pdgcode==-13 || pdgcode==13) ptmu = pt;             
156         }
157
158 // global trigger
159         TString tree("D");
160         if ( readFromRP ) tree = "R";
161         
162         AliMUONVTriggerStore* triggerStore = diRec.TriggerStore(ievent,tree.Data());
163         AliMUONGlobalTrigger* gloTrg = triggerStore->Global();
164         
165         if (gloTrg->SingleLpt()>=1 ) {
166             lptmuon++;
167             ptlpt->Fill(ptmu);
168         }       
169         if (gloTrg->SingleHpt()>=1 ) {
170             hptmuon++;
171             pthpt->Fill(ptmu);
172         }
173
174 // Hits
175         Int_t itrack, ntracks, NbHits[4];
176         Int_t SumNbHits;
177         
178         for (Int_t j=0; j<4; j++) NbHits[j]=0;
179  
180         ntracks = (Int_t) diSim.NumberOfTracks(ievent);
181    
182         for (itrack=0; itrack<ntracks; itrack++) { // track loop
183             AliMUONVHitStore* hitStore = diSim.HitStore(ievent,itrack);      
184             AliMUONHit* mHit;
185             TIter next(hitStore->CreateIterator());
186
187             while ( ( mHit = static_cast<AliMUONHit*>(next()) ) )
188             {
189                 Int_t Nch        = mHit->Chamber(); 
190                 Int_t hittrack   = mHit->Track();
191                 Float_t IdPart   = mHit->Particle();        
192             
193                 for (Int_t j=0;j<4;j++) {
194                     Int_t kch=11+j;
195                     if (hittrack==0) {
196                         if (Nch==kch && (IdPart==-13 || IdPart==13) && NbHits[j]==0) NbHits[j]++; 
197                     }               
198                 }    
199             }
200         } // end track loop
201
202
203 // 3/4 coincidence 
204         SumNbHits=NbHits[0]+NbHits[1]+NbHits[2]+NbHits[3];
205
206         if (SumNbHits==3 || SumNbHits==4) {
207             coincmuon++;
208             ptcoinc->Fill(ptmu);
209         }            
210                 
211
212     } // end loop on event
213
214     if (coincmuon==0) {
215         cout << " >>> <E> coincmuon = 0 after event loop " << "\n";
216         cout << " >>> this probably means that input does not contain one (and only one) " << "\n";
217         cout << " >>> muon track per event as it should " << "\n";
218         cout << " >>> see README for further information " << "\n";
219         cout << " >>> exiting now ! " << "\n";
220         return;
221     }
222
223
224     fprintf(fdat,"\n");    
225     fprintf(fdat,"\n");
226     fprintf(fdat," Number of events = %d \n",nevents);  
227     fprintf(fdat," Number of events with 3/4 coinc = %d \n",(Int_t)coincmuon);
228     fprintf(fdat," Nomber of dimuons with 3/4 coinc Lpt cut = %d \n",(Int_t)lptmuon); 
229     fprintf(fdat," Number of dimuons with 3/4 coinc Hpt cut = %d \n",(Int_t)hptmuon);  
230     fprintf(fdat,"\n"); 
231     
232     Double_t efficiency,error;
233     
234     efficiency=lptmuon/coincmuon;  
235     error=efficiency*TMath::Sqrt((lptmuon+coincmuon)/(lptmuon*coincmuon));  
236     fprintf(fdat," Efficiency Lpt cut = %4.4f +/- %4.4f\n",efficiency,error);
237     
238     efficiency=hptmuon/coincmuon; 
239     error=efficiency*TMath::Sqrt((hptmuon+coincmuon)/(hptmuon*coincmuon)); 
240     fprintf(fdat," Efficiency Hpt cut = %4.4f +/- %4.4f\n",efficiency,error);
241     
242     fclose(fdat);
243     
244     Double_t x1,x2,xval,xerr;
245     
246     for (Int_t i=0;i<50;i++) {
247         x1=ptcoinc->GetBinContent(i+1);
248         
249         x2=ptlpt->GetBinContent(i+1);
250         if (x1!=0 && x2!=0) {
251             xval=x2/x1;
252             xerr=xval*TMath::Sqrt((x1+x2)/x1*x2);   
253             ptlpt->SetBinContent(i+1,xval);
254             ptlpt->SetBinError(i+1,0);
255         } else {
256             ptlpt->SetBinContent(i+1,0.);
257             ptlpt->SetBinError(i+1,0.);     
258         } 
259         
260         x2=pthpt->GetBinContent(i+1);
261         if (x1!=0 && x2!=0) {
262             xval=x2/x1;
263             xerr=xval*TMath::Sqrt((x1+x2)/x1*x2);   
264             pthpt->SetBinContent(i+1,xval);
265             pthpt->SetBinError(i+1,0);
266         } else {
267             pthpt->SetBinContent(i+1,0.);
268             pthpt->SetBinError(i+1,0.);     
269         }              
270     }      
271     
272     TF1 *fitlpt = new TF1("fitlpt",fitArc,0.,5.,4); 
273     TF1 *fithpt = new TF1("fithpt",fitArc,0.,5.,4);      
274     
275     TCanvas *c1 = new TCanvas("c1","Trigger efficiency vs pt and y for single muon",200,0,900,400); 
276     c1->Divide(2,1);      
277     
278     c1->cd(1);
279     ptlpt->SetMinimum(0.);
280     ptlpt->SetMaximum(1.05);   
281     ptlpt->SetTitle("");
282     ptlpt->GetXaxis()->SetTitle("P_{T} (GeV/c)");
283     ptlpt->GetYaxis()->SetTitle("Efficiency");     
284     //ptlpt->GetXaxis()->SetRange(3);
285     ptlpt->Draw(""); 
286     fitlpt->SetParameters(0.602,0.774,0.273,0.048);  
287     fitlpt->SetLineColor(2);
288     fitlpt->SetLineWidth(2);
289     fitlpt->Draw("SAME");    
290     TLegend * leg = new TLegend(0.5,0.38,0.85,0.53); 
291     leg = new TLegend(0.5,0.38,0.85,0.53); 
292     leg->SetBorderSize(0);
293     leg->SetFillColor(0);
294     leg->SetTextSize(0.05);
295     leg->SetTextFont(22);
296     leg->SetHeader("Lpt trigger pt cut");
297     leg->AddEntry(fitlpt,"Ref.","l");
298     leg->AddEntry(ptlpt,"New","l");                 
299     leg->Draw("SAME");
300     
301     c1->cd(2);
302     pthpt->SetMinimum(0.);
303     pthpt->SetMaximum(1.05);   
304     pthpt->SetTitle("");
305     pthpt->GetXaxis()->SetTitle("P_{T} (GeV/c)");  
306     pthpt->GetYaxis()->SetTitle("Efficiency");     
307     //pthpt->GetXaxis()->SetRange(3);
308     pthpt->Draw(""); 
309     fithpt->SetParameters(0.627,0.393,0.855,0.0081); 
310     fithpt->SetLineColor(2);
311     fithpt->SetLineWidth(2);
312     fithpt->Draw("SAME");    
313     leg = new TLegend(0.5,0.38,0.85,0.53); 
314     leg->SetBorderSize(0);
315     leg->SetFillColor(0);
316     leg->SetTextSize(0.05);
317     leg->SetTextFont(22);
318     leg->SetHeader("Hpt trigger pt cut");
319     leg->AddEntry(fithpt,"Ref.","l");
320     leg->AddEntry(pthpt,"New","l");                 
321     leg->Draw("SAME");
322     
323     c1->SaveAs("MUONTriggerEfficiencyPt.gif");
324     c1->SaveAs("MUONTriggerEfficiencyPt.eps");      
325
326 }