]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/muon/AliAnalysisTaskSingleMuESD.cxx
Updates to run with deltas (L. Cunqueiro)
[u/mrichter/AliRoot.git] / PWG / muon / AliAnalysisTaskSingleMuESD.cxx
1
2 /* $Id$ */
3
4 #include "TChain.h"
5 #include "TTree.h"
6 #include "TNtuple.h"
7 #include "TH1F.h"
8 #include "TString.h"
9
10 #include "AliAnalysisTask.h"
11 #include "AliAnalysisManager.h"
12
13 #include "AliESDEvent.h"
14 #include "AliESDInputHandler.h"
15 #include "AliESDVertex.h"
16 #include "AliESDMuonTrack.h"
17
18 #include "AliAnalysisTaskSingleMuESD.h"
19
20 // analysis task for single muon analysis of the ESD events
21 // Authors: Bogdan Vulpescu, Nicole Bastid, LPC Clermont-Ferrand
22
23 ClassImp(AliAnalysisTaskSingleMuESD)
24
25 //________________________________________________________________________
26 AliAnalysisTaskSingleMuESD::AliAnalysisTaskSingleMuESD(const char *name) 
27   : AliAnalysisTask(name, ""), 
28     fESD(0), 
29     fNtuple(0), 
30     fTrigger(0), 
31     fMaskTrig1MuL(1),
32     fMaskTrig1MuH(2),
33     fMaskTrigUSL(4),
34     fMaskTrigLSL(16),
35     fMaskTrigUSH(8),
36     fMaskTrigLSH(32),
37     fTriggerType(kFALSE)
38 {
39   // Constructor
40
41   // Define input and output slots here
42   // Input slot #0 works with a TChain
43   DefineInput(0, TChain::Class());
44   // Output slot #0 writes into a TNtuple container
45   DefineOutput(0, TNtuple::Class());
46   // Output slot #1 writes into a TH1F container
47   DefineOutput(1, TH1F::Class());
48 }
49
50 //________________________________________________________________________
51 void AliAnalysisTaskSingleMuESD::ConnectInputData(Option_t *) 
52 {
53   // Connect ESD or AOD here
54   // Called once
55
56   TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
57   if (!tree) {
58     Printf("ERROR: Could not read chain from input slot 0");
59   } else {
60     // Disable all branches and enable only the needed ones
61     // The next two lines are different when data produced as AliESDEvent is read
62     tree->SetBranchStatus("*", kFALSE);
63     tree->SetBranchStatus("MuonTracks.*", kTRUE);
64     tree->SetBranchStatus("SPDVertex.*", kTRUE);
65     tree->SetBranchStatus("AliESDHeader.*", kTRUE);
66
67     AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
68
69     if (!esdH) {
70       Printf("ERROR: Could not get ESDInputHandler");
71     } else
72       fESD = esdH->GetEvent();
73   }
74 }
75
76 //________________________________________________________________________
77 void AliAnalysisTaskSingleMuESD::CreateOutputObjects()
78 {
79   // Create output ntuple
80
81   fNtuple = new TNtuple("ntEsd1mu","ntEsd1mu","VtZ:Mult:TMask:Z:Pt:P:Eta:Theta:Phi:Zpos:Hits:Xi2:Eff:TrigMatch:TrigXi2");
82   // histogram with trigger classes
83   fTrigger = new TH1F("TriggerMask","L/H/USL/USH/LSL/LSH",6,0,6);
84
85 }
86
87 //________________________________________________________________________
88 void AliAnalysisTaskSingleMuESD::Exec(Option_t *) 
89 {
90   // Main loop
91   // Called for each event
92   // Fill the ntuple with variables and the trigger histogram
93
94   Float_t  ntvar[15], pt, p, eta, thetad, phid, charge, chi2, zcoor, chi2MatchTrig;
95   ULong64_t mask = 0;
96   Int_t masksave = 0, hits, matchTrig;
97   Double_t rEff[2];
98   Double_t fZVertex = 0;
99   Double_t fYVertex = 0;
100   Double_t fXVertex = 0;
101   Double_t fErXVertex = 0;
102   Double_t fErYVertex = 0;
103   
104   if (!fESD) {
105     Printf("ERROR: fESD not available");
106     return;
107   }
108
109   // Trigger checks
110   // mask values: L/H/USL/USH/LSL/LSH
111   // Bit position    Histogram bin     Trigger class
112   // 1 (less)        1 ("L")           MULow          single low pt
113   // 2               2 ("H")           MUHigh         single high pt
114   // 3               3 ("USL")         MULU           unlike-sign low pt
115   // 4               4 ("USH")         MUHU           unlike-sign high pt
116   // 5               5 ("LSL")         MULL           like-sign low pt
117   // 6 (most)        6 ("LSH")         MUHL           like-sign high pt
118   if (fTriggerType) {
119     mask = fESD->GetTriggerMask();
120     if (mask & fMaskTrig1MuL) {
121       fTrigger->Fill(0.5);
122       masksave |= 0x01;
123     }
124     if (mask & fMaskTrig1MuH) {
125       fTrigger->Fill(1.5);
126       masksave |= 0x02;
127     }
128     if (mask & fMaskTrigUSL)  {
129       fTrigger->Fill(2.5);
130       masksave |= 0x04;
131     }
132     if (mask & fMaskTrigUSH)  {
133       fTrigger->Fill(3.5);
134       masksave |= 0x08;
135     }
136     if (mask & fMaskTrigLSL)  {
137       fTrigger->Fill(4.5);
138       masksave |= 0x10;
139     }
140     if (mask & fMaskTrigLSH)  {
141       fTrigger->Fill(5.5);
142       masksave |= 0x20;
143     }
144   } else {
145     TString firedTC = fESD->GetFiredTriggerClasses();
146     if (firedTC.Contains("MULow")) {
147       fTrigger->Fill(0.5);
148       masksave |= 0x01;
149     }
150     if (firedTC.Contains("MUHigh")) {
151       fTrigger->Fill(1.5);
152       masksave |= 0x02;
153     }
154     if (firedTC.Contains("MULU")) {
155       fTrigger->Fill(2.5);
156       masksave |= 0x04;
157     }
158     if (firedTC.Contains("MUHU")) {
159       fTrigger->Fill(3.5);
160       masksave |= 0x08;
161     }
162     if (firedTC.Contains("MULL")) {
163       fTrigger->Fill(4.5);
164       masksave |= 0x10;
165     }
166     if (firedTC.Contains("MUHL")) {
167       fTrigger->Fill(5.5);
168       masksave |= 0x20;
169     }
170   }
171
172   // get the SPD reconstructed vertex (vertexer) 
173   AliESDVertex* vertex = (AliESDVertex*) fESD->GetVertex();
174   if (vertex->GetNContributors()) {
175     fZVertex = vertex->GetZv();
176     fYVertex = vertex->GetYv();
177     fXVertex = vertex->GetXv();
178     fErXVertex = vertex->GetXRes();
179     fErYVertex = vertex->GetYRes();
180   }
181
182   Int_t nTracks = (Int_t)fESD->GetNumberOfMuonTracks() ;
183   for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
184
185     AliESDMuonTrack* muonTrack = new AliESDMuonTrack(*(fESD->GetMuonTrack(iTrack)));
186     if (!muonTrack) {
187       Printf("ERROR: Could not receive track %d", iTrack);
188       continue;
189     }
190
191     pt = muonTrack->Pt();
192     p = muonTrack->P();
193     eta = muonTrack->Eta();
194     thetad =  muonTrack->Theta()*180./TMath::Pi();      
195     phid = muonTrack->Phi()*180./TMath::Pi();
196     charge = muonTrack->Charge();
197     hits = muonTrack->GetNHit();        
198     chi2 = muonTrack->GetChi2() / (2.0 * hits - 5);
199     zcoor = muonTrack->GetZ();
200     matchTrig = muonTrack->GetMatchTrigger();
201     chi2MatchTrig = muonTrack->GetChi2MatchTrigger();
202     
203     rEff[0] = rEff[1] = 0.0;
204     GetEffFitted(eta,pt,rEff);
205         
206     ntvar[0] = fZVertex;
207     ntvar[1] = nTracks;
208     ntvar[2] = masksave;
209     ntvar[3] = charge;
210     ntvar[4] = pt;
211     ntvar[5] = p;
212     ntvar[6] = eta;     
213     ntvar[7] = thetad;
214     ntvar[8] = phid;
215     ntvar[9] = zcoor;
216     ntvar[10] = hits;   
217     ntvar[11] = chi2;
218     ntvar[12] = rEff[0];
219     ntvar[13] = matchTrig;
220     ntvar[14] = chi2MatchTrig;  
221
222     fNtuple->Fill(ntvar);
223
224
225   } //track loop 
226
227   // Post output data.
228   PostData(0, fNtuple);
229   PostData(1, fTrigger);
230
231 }      
232
233 //________________________________________________________________________
234 void AliAnalysisTaskSingleMuESD::Terminate(Option_t *) 
235 {
236   // the end
237 }
238
239 //________________________________________________________________________
240 void AliAnalysisTaskSingleMuESD::GetEffFitted(Double_t eta, Double_t pt, Double_t rEff[2])
241 {
242   // tracking efficiency
243
244   Int_t nbinEta =15;
245   Float_t etaMin = -4.;
246   Float_t etaBin = 0.1;
247   Int_t numBin = -1;
248   rEff[0]=1.;rEff[1]=0.;
249   
250   // nbinEta = number of eta bins startig at etaMin = -4 and binning of 0.1 (etaBin)
251   Float_t p0[15] = {-0.9338,-0.9441,-0.9397,-0.9586,-0.9713,-0.9752,-0.9711,-0.9751,-0.9746,-0.9729,-0.9754,-0.9722,-0.9722,-0.971,-0.9709};
252   Float_t p1[15] = {-0.6383,-1.4729,-1.0109,-1.4316,-1.5926,-1.367,-1.1895,-1.2834,-1.3289,-1.5916,-1.4258,-1.0983,-1.0812,-0.7179, -0.4613};      
253   
254   Float_t etaInf, etaSup;
255   for (Int_t i = 1; i < nbinEta+1; i++){
256     etaInf = etaMin + (i-1)*etaBin;
257     etaSup = etaMin + i*etaBin;
258     if(eta > etaInf && eta <= etaSup) numBin = i-1;
259   } 
260   
261   if(numBin>=0 && pt >=1.) rEff[0] = p0[numBin] * TMath::TanH(p1[numBin]*pt);
262   rEff[1] = 0;
263
264 }
265
266 //________________________________________________________________________
267 void AliAnalysisTaskSingleMuESD::SetTriggerType(const Char_t *trig)
268 {
269   // set the trigger masks according to the trigger type used in the
270   // simulations
271   
272   TString triggerType(trig);
273
274   if (!triggerType.CompareTo("MUON")) {
275     fMaskTrig1MuL =  1;  // MULow
276     fMaskTrig1MuH =  2;  // MUHigh
277     fMaskTrigUSL  =  4;  // MULU
278     fMaskTrigLSL  = 16;  // MULL
279     fMaskTrigUSH  =  8;  // MUHU
280     fMaskTrigLSH  = 32;  // MUHL
281   } else if (!triggerType.CompareTo("p-p")) {
282     fMaskTrig1MuL =    128;  // MULow
283     fMaskTrig1MuH =    256;  // MUHigh
284     fMaskTrigUSL  =  65536;  // MULU
285     fMaskTrigLSL  = 131072;  // MULL
286     fMaskTrigUSH  = 262144;  // MUHU
287     fMaskTrigLSH  = 524288;  // MUHL
288   } else {
289     // MUON trigger, values by default in the constructor
290   }
291   
292   fTriggerType = kTRUE;
293
294 }
295