]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONTriggerEfficiencyPt.C
Adding comments only
[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    3 pt cuts. Results are compared to the reference (red curves).   
21    To be used with (at least) 10000 events as follows
22    AliGenBox * gener = new AliGenBox(1);
23    gener->SetPtRange(0.,10.);
24    gener->SetPhiRange(0., 360.);         
25    gener->SetThetaRange(171.000,178.001);
26    gener->SetPart(13);           // or -13
27    gener->SetOrigin(0.,0., 0.);  
28    gener->SetSigma(0.0, 0.0, 0.0);     
29
30    Author: Fabien Guerin, LPC Clermont-Ferrand, Jan. 2006
31 */
32
33
34 // ROOT includes
35 #include "TBranch.h"
36 #include "TClonesArray.h"
37 #include "TFile.h"
38 #include "TF1.h"
39 #include "TH1.h"
40 #include "TH2.h"
41 #include "TParticle.h"
42 #include "TTree.h"
43 #include "TMath.h"
44 #include "TStyle.h"
45 #include "TCanvas.h"
46 #include "TLegend.h"
47
48
49 // STEER includes
50 #include "AliRun.h"
51 #include "AliDetector.h"
52 #include "AliRunLoader.h"
53 #include "AliHeader.h"
54 #include "AliLoader.h"
55 #include "AliStack.h"
56 #include "AliSegmentation.h"
57 #include "AliMC.h"
58
59 // MUON includes
60 #include "AliMUON.h"
61 #include "AliMUONData.h"
62 #include "AliMUONHit.h"
63 #include "AliMUONChamber.h"
64 #include "AliMUONConstants.h"
65 #include "AliMUONDigit.h"
66 #include "AliMUONRawCluster.h"
67 #include "AliMUONGlobalTrigger.h"
68 #include "AliMUONLocalTrigger.h"
69 #include "AliMUONTrack.h"
70
71 Double_t fitArc(Double_t *x,Double_t *par) 
72 {
73   Double_t h1=par[1]*(x[0]+par[2]);
74   return par[0]*TMath::ATan(TMath::Power(h1,6))+par[3];
75 }
76
77 Double_t fitArch(Double_t *x,Double_t *par) 
78 {
79   Double_t h0=2*par[0]/(TMath::Pi());
80   Double_t h1=par[1]*x[0]+par[2]*x[0]*x[0];
81   return h0*TMath::TanH(h1)+par[3];
82 }
83
84 void MUONTriggerEfficiencyPt(char filename[10]="galice.root")
85 {
86
87 // define style
88     TStyle *st1 = new TStyle("st1","My STYLE");
89     
90     st1->SetOptStat(0);
91     st1->SetOptFit(111);  
92     st1->SetFrameFillColor(10);
93     st1->SetCanvasColor(10);
94     st1->SetFillColor(10);
95     st1->SetTitleBorderSize(0);
96     st1->SetTitleOffset(1.2,"XY");
97     st1->SetTitleSize(0.05,"XY"); 
98     st1->SetLabelSize(0.045,"XY");  
99     st1->SetLabelFont(22,"XY");  
100     st1->SetTitleFont(22,"XY");     
101     st1->SetOptTitle(0);   
102     st1->SetStatColor(10);
103     st1->SetStatFont(62);     
104     st1->SetFitFormat("4.4g");
105     st1->SetPadLeftMargin(0.15);
106     st1->SetPadRightMargin(0.06); 
107     st1->SetPadTopMargin(0.06);
108     st1->SetPadBottomMargin(0.15); 
109     st1->cd();
110     
111     gROOT->ForceStyle();
112     //TGaxis::SetMaxDigits(3);
113     
114 // beginning of macro    
115     TH1F *ptlpt    = new TH1F("ptlpt","",50,0.15,5.); 
116     TH1F *pthpt    = new TH1F("pthpt","",50,0.15,5.);
117     TH1F *ptcoinc  = new TH1F("ptcoinc","",50,0.15,5.);  
118         
119     TParticle *particle;
120     AliStack* stack; 
121     
122     Double_t coincmuon=0;
123     Double_t lptmuon=0;
124     Double_t hptmuon=0;
125     Double_t ptmu=0;
126
127 // output file
128     char digitdat[100];
129     sprintf(digitdat,"MUONTriggerEfficiencyPt.out");  
130     FILE *fdat=fopen(digitdat,"w");
131
132 // Creating Run Loader and openning file containing Hits
133     AliRunLoader * RunLoader = AliRunLoader::Open(filename,"MUONFolder","READ");
134     
135     if (RunLoader ==0x0) {
136         printf(">>> Error : Error Opening %s file \n",filename);
137         return;
138     }
139     
140     Int_t nevents = RunLoader->GetNumberOfEvents();
141     
142 // loading hits
143     AliLoader * MUONLoader = RunLoader->GetLoader("MUONLoader");     
144     MUONLoader->LoadHits("READ");
145     AliMUONData data_hits(MUONLoader,"MUON","MUON");
146     
147     TClonesArray * globalTrigger;
148     AliMUONGlobalTrigger * gloTrg;
149
150 // loading trigger
151     MUONLoader->LoadDigits("READ");  
152     AliMUONData muondata(MUONLoader,"MUON","MUON");
153
154 // Loading kine
155     RunLoader->LoadKinematics("READ");
156     
157     for (Int_t ievent=0; ievent<nevents; ievent++) {  // Event loop
158         
159         RunLoader->GetEvent(ievent);
160         
161         if (ievent%500==0) printf("ievent = %d \n",ievent);
162
163
164 // kine
165         Int_t iparticle, nparticles;
166         stack = RunLoader->Stack();
167         nparticles = (Int_t) stack->GetNtrack();        
168         for (iparticle=0; iparticle<nparticles; iparticle++) {             
169             particle = stack->Particle(iparticle);   
170             Float_t pt = particle->Pt();            
171             Int_t pdgcode = particle->GetPdgCode();                        
172             if (pdgcode==-13 || pdgcode==13) ptmu = pt;             
173         }
174
175 // trigger 
176         muondata.SetTreeAddress("D,GLT");
177         muondata.GetTriggerD();
178     
179         globalTrigger = muondata.GlobalTrigger();
180
181         Int_t nglobals = (Int_t) globalTrigger->GetEntriesFast(); // should be 1 
182         for (Int_t iglobal=0; iglobal<nglobals; iglobal++) { // Global Trigger
183             gloTrg = static_cast<AliMUONGlobalTrigger*>(globalTrigger->At(iglobal));
184             
185             if (gloTrg->SingleLpt()>=1 ) {
186                 lptmuon++;
187                 ptlpt->Fill(ptmu);
188             }       
189             if (gloTrg->SingleHpt()>=1 ) {
190                 hptmuon++;
191                 pthpt->Fill(ptmu);
192             }
193         } // end of loop on Global Trigger      
194         muondata.ResetTrigger();  
195
196 // Hits
197         RunLoader->GetEvent(ievent);  
198         data_hits.SetTreeAddress("H");    
199         
200         Int_t itrack, ntracks, NbHits[4];
201         Int_t SumNbHits;
202
203         for (Int_t j=0; j<4; j++) NbHits[j]=0;
204  
205         ntracks = (Int_t) data_hits.GetNtracks();      
206      
207         for (itrack=0; itrack<ntracks; itrack++) { // track loop
208          data_hits.GetTrack(itrack); 
209
210          Int_t ihit, nhits;
211          nhits = (Int_t) data_hits.Hits()->GetEntriesFast();   
212          AliMUONHit* mHit;
213
214     
215           for (ihit=0; ihit<nhits; ihit++) {
216             mHit = static_cast<AliMUONHit*>(data_hits.Hits()->At(ihit));
217             Int_t Nch        = mHit->Chamber(); 
218             Int_t hittrack   = mHit->Track();
219             Int_t IdPart     = mHit->Particle();            
220
221             for (Int_t j=0;j<4;j++) {
222               Int_t kch=11+j;
223               if (hittrack==0) {
224                 if (Nch==kch && (IdPart==-13 || IdPart==13) && NbHits[j]==0) NbHits[j]++; 
225                 }               
226             }    
227           }
228           data_hits.ResetHits();
229         } // end track loop
230         
231 // 3/4 coincidence 
232         SumNbHits=NbHits[0]+NbHits[1]+NbHits[2]+NbHits[3];
233
234         if (SumNbHits==3 || SumNbHits==4) {
235             coincmuon++;
236             ptcoinc->Fill(ptmu);
237         }            
238                 
239       } // end loop on event
240       
241       MUONLoader->UnloadHits();
242       MUONLoader->UnloadDigits();
243       MUONLoader->UnloadRecPoints();
244       RunLoader->UnloadKinematics();   
245
246       delete RunLoader;
247      
248       fprintf(fdat,"\n");    
249       fprintf(fdat,"\n");
250       fprintf(fdat," Number of events = %d \n",nevents);  
251       fprintf(fdat," Number of events with 3/4 coinc = %d \n",(Int_t)coincmuon);
252       fprintf(fdat," Nomber of dimuons with 3/4 coinc Lpt cut = %d \n",(Int_t)lptmuon); 
253       fprintf(fdat," Number of dimuons with 3/4 coinc Hpt cut = %d \n",(Int_t)hptmuon);  
254       fprintf(fdat,"\n"); 
255       
256       Double_t efficiency,error;
257
258       efficiency=lptmuon/coincmuon;  
259       error=efficiency*TMath::Sqrt((lptmuon+coincmuon)/(lptmuon*coincmuon));  
260       fprintf(fdat," Efficiency Lpt cut = %4.4f +/- %4.4f\n",efficiency,error);
261       
262       efficiency=hptmuon/coincmuon; 
263       error=efficiency*TMath::Sqrt((hptmuon+coincmuon)/(hptmuon*coincmuon)); 
264       fprintf(fdat," Efficiency Hpt cut = %4.4f +/- %4.4f\n",efficiency,error);
265       
266       fclose(fdat);
267       
268       Double_t x1,x2,xval,xerr;
269       
270       for (Int_t i=0;i<50;i++) {
271           x1=ptcoinc->GetBinContent(i+1);
272           
273           x2=ptlpt->GetBinContent(i+1);
274           if (x1!=0 && x2!=0) {
275               xval=x2/x1;
276               xerr=xval*TMath::Sqrt((x1+x2)/x1*x2);   
277               ptlpt->SetBinContent(i+1,xval);
278               ptlpt->SetBinError(i+1,0);
279           } else {
280               ptlpt->SetBinContent(i+1,0.);
281               ptlpt->SetBinError(i+1,0.);     
282           } 
283           
284           x2=pthpt->GetBinContent(i+1);
285           if (x1!=0 && x2!=0) {
286               xval=x2/x1;
287               xerr=xval*TMath::Sqrt((x1+x2)/x1*x2);   
288               pthpt->SetBinContent(i+1,xval);
289               pthpt->SetBinError(i+1,0);
290           } else {
291               pthpt->SetBinContent(i+1,0.);
292               pthpt->SetBinError(i+1,0.);     
293           }              
294       }
295       
296       
297       TF1 *fitlpt = new TF1("fitlpt",fitArc,0.,5.,4); 
298       TF1 *fithpt = new TF1("fithpt",fitArc,0.,5.,4);      
299       
300       TCanvas *c1 = new TCanvas("c1","Trigger efficiency vs pt and y for single muon",200,0,900,400); 
301       c1->Divide(2,1);      
302       
303       c1->cd(1);
304       ptlpt->SetMinimum(0.);
305       ptlpt->SetMaximum(1.05);   
306       ptlpt->SetTitle("");
307       ptlpt->GetXaxis()->SetTitle("P_{T} (GeV/c)");
308       ptlpt->GetYaxis()->SetTitle("Efficiency");     
309       //ptlpt->GetXaxis()->SetRange(3);
310       ptlpt->Draw(""); 
311       fitlpt->SetParameters(0.602,0.774,0.273,0.048);  
312       fitlpt->SetLineColor(2);
313       fitlpt->SetLineWidth(2);
314       fitlpt->Draw("SAME");    
315       TLegend * leg = new TLegend(0.5,0.38,0.85,0.53); 
316       leg = new TLegend(0.5,0.38,0.85,0.53); 
317       leg->SetBorderSize(0);
318       leg->SetFillColor(0);
319       leg->SetTextSize(0.05);
320       leg->SetTextFont(22);
321       leg->SetHeader("Lpt trigger pt cut");
322       leg->AddEntry(fitlpt,"Ref.","l");
323       leg->AddEntry(ptlpt,"New","l");                 
324       leg->Draw("SAME");
325       
326       c1->cd(2);
327       pthpt->SetMinimum(0.);
328       pthpt->SetMaximum(1.05);   
329       pthpt->SetTitle("");
330       pthpt->GetXaxis()->SetTitle("P_{T} (GeV/c)");  
331       pthpt->GetYaxis()->SetTitle("Efficiency");     
332       //pthpt->GetXaxis()->SetRange(3);
333       pthpt->Draw(""); 
334       fithpt->SetParameters(0.627,0.393,0.855,0.0081); 
335       fithpt->SetLineColor(2);
336       fithpt->SetLineWidth(2);
337       fithpt->Draw("SAME");    
338       leg = new TLegend(0.5,0.38,0.85,0.53); 
339       leg->SetBorderSize(0);
340       leg->SetFillColor(0);
341       leg->SetTextSize(0.05);
342       leg->SetTextFont(22);
343       leg->SetHeader("Hpt trigger pt cut");
344       leg->AddEntry(fithpt,"Ref.","l");
345       leg->AddEntry(pthpt,"New","l");                 
346       leg->Draw("SAME");
347       
348       c1->SaveAs("MUONTriggerEfficiencyPt.gif");
349       c1->SaveAs("MUONTriggerEfficiencyPt.eps");
350       
351 }