]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/muon/AliAnalysisTaskMuonCollisionMultiplicity.cxx
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PWG / muon / AliAnalysisTaskMuonCollisionMultiplicity.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 /// Compute the number of Muons tracks as a function of the SPD tracklets multiplicity
17 /// Compare with Monte Carlo tracks
18 /// Author Matthieu LENHARDT - SUBATECH, Nantes
19
20
21 //PWG3/muon includes
22 #include "AliAnalysisTaskMuonCollisionMultiplicity.h"
23
24 //STEER includes
25 #include "AliAODEvent.h"
26 #include "AliAODVertex.h"
27 #include "AliESDEvent.h"
28 #include "AliESDVertex.h"
29 #include "AliAODTrack.h"
30 #include "AliESDMuonTrack.h"
31 #include "AliAODTracklets.h"
32 #include "AliAODInputHandler.h"
33 #include "AliAnalysisManager.h"
34 #include "AliAODDimuon.h"
35 #include "AliMCEvent.h"
36 #include "AliMCParticle.h"
37 #include "AliMultiplicity.h"
38
39 //ROOT includes
40 #include <TH2D.h>
41 #include <THnSparse.h>
42 #include <TChain.h>
43 #include <TList.h>
44 #include <TArrayD.h>
45 #include <Riostream.h>
46 #include <TParticle.h>
47 #include <TLorentzVector.h>
48
49 //___________________________________________________
50 AliAnalysisTaskMuonCollisionMultiplicity::AliAnalysisTaskMuonCollisionMultiplicity()
51   :
52   AliAnalysisTaskSE(),
53   fIsInit(kFALSE),
54   fAOD(0x0),
55   fESD(0x0),
56   fEtaCut(0),
57   fTrackletMultiplicity(0),
58   fTriggerList(0),
59   fSingleMuonList(0),
60   fDimuonList(0),
61   fMonteCarloList(0)
62 {
63   ///Default Constructor
64 }
65
66
67 //___________________________________________________
68 AliAnalysisTaskMuonCollisionMultiplicity::AliAnalysisTaskMuonCollisionMultiplicity(const AliAnalysisTaskMuonCollisionMultiplicity& src)
69   :
70   AliAnalysisTaskSE(),
71   fIsInit(kFALSE),
72   fAOD(0x0),
73   fESD(0x0),
74   fEtaCut(0),
75   fTrackletMultiplicity(0),
76   fTriggerList(0),
77   fSingleMuonList(0),
78   fDimuonList(0),
79   fMonteCarloList(0)
80 {
81   /// copy ctor
82   src.Copy(*this);
83
84 }
85
86 //___________________________________________________
87 AliAnalysisTaskMuonCollisionMultiplicity& AliAnalysisTaskMuonCollisionMultiplicity::operator=(const AliAnalysisTaskMuonCollisionMultiplicity& src)
88 {
89   /// assignement operator
90   if ( this != &src ) 
91   {
92     src.Copy(*this);
93   }
94   return *this;
95 }
96
97
98
99
100 //___________________________________________________
101 AliAnalysisTaskMuonCollisionMultiplicity::AliAnalysisTaskMuonCollisionMultiplicity(const Char_t *name)
102   :
103   AliAnalysisTaskSE(name),
104   fIsInit(kFALSE),
105   fAOD(0x0),
106   fESD(0x0),
107   fEtaCut(0),
108   fTrackletMultiplicity(0),
109   fTriggerList(0),
110   fSingleMuonList(0),
111   fDimuonList(0),
112   fMonteCarloList(0)
113 {
114   // Define Inputs and outputs
115   //DefineInput(0, TChain::Class());
116   //DefineInput(1, TChain::Class());
117   DefineOutput(0, TList::Class());
118   DefineOutput(1, TList::Class());
119   DefineOutput(2, TList::Class());
120   DefineOutput(3, TList::Class());
121   DefineOutput(4, TList::Class());
122 }
123
124
125 //______________________________________________________________________________
126 AliAnalysisTaskMuonCollisionMultiplicity::~AliAnalysisTaskMuonCollisionMultiplicity()
127 {
128 // Destructor.
129   delete fAOD;
130   delete fESD;
131   delete fTriggerList;
132   delete fSingleMuonList;
133   delete fDimuonList;
134   delete fMonteCarloList;
135 }
136
137
138
139 //________________________________________________________________________
140 void AliAnalysisTaskMuonCollisionMultiplicity::UserCreateOutputObjects()
141 {
142   // Initialise the object, and open the file.
143   if (!fIsInit)
144     Init();
145   OpenFile(0);
146 }
147
148
149
150
151 //________________________________________________________________________
152 void AliAnalysisTaskMuonCollisionMultiplicity::UserExec(Option_t */*option*/)
153 {
154   // Execute the analysis task
155   fAOD = 0x0;
156   fESD = 0x0;
157
158   if (!fIsInit)
159     Init();
160
161   fAOD = dynamic_cast<AliAODEvent *> (InputEvent());
162   if (!fAOD)
163     fESD = dynamic_cast<AliESDEvent *> (InputEvent());
164   
165  
166   if (fAOD)
167     CheckEventAOD();
168
169   if (fESD)
170     CheckEventESD();
171
172   PostData(1, fTriggerList);
173   PostData(2, fSingleMuonList);
174   PostData(3, fDimuonList);
175   PostData(4, fMonteCarloList);
176 }
177
178
179
180 //________________________________________________________________________
181 void AliAnalysisTaskMuonCollisionMultiplicity::NotifyRun()
182 {
183   // Notify run
184 }
185
186
187 //________________________________________________________________________
188 void AliAnalysisTaskMuonCollisionMultiplicity::FinishTaskOutput()
189 {
190   // Finish the task
191 }
192
193
194 //__________________________________________________________________________
195 Bool_t AliAnalysisTaskMuonCollisionMultiplicity::CheckEventAOD()
196 {
197   // Check if the AOD event pass the cuts
198   AliAODVertex *vertex = fAOD->GetPrimaryVertex();
199   
200   if (!vertex)
201      return kFALSE;
202
203   if (vertex->GetNContributors() < 1)
204     return kFALSE;
205   
206   ComputeMultiplicity();
207   
208   // Variables use to determine the type of trigger :
209   // 0 for minimum bias : CINT1B, CINT1-B or MB1
210   // 1 for muon events : CMUS1B, CMUS1-B or MULow
211   // -1 for everything else
212   TString trigger = fAOD->GetFiredTriggerClasses();
213   Int_t triggerClass = -1; 
214
215   if (trigger.Contains("CINT1B") || trigger.Contains("CINT1-B") || trigger.Contains("MB1"))
216     triggerClass = 0;
217
218   if (trigger.Contains("CMUS1B") || trigger.Contains("CMUS1-B") || trigger.Contains("MULow"))
219     triggerClass = 1;
220
221   if (triggerClass >= 0)
222     FillHistosAOD(triggerClass);
223   
224   return kTRUE;
225 }
226
227
228
229 //__________________________________________________________________________
230 Bool_t AliAnalysisTaskMuonCollisionMultiplicity::CheckEventESD()
231 {
232   // Check if the ESD event pass the cuts
233   const AliESDVertex *vertex = fESD->GetPrimaryVertex();
234   
235   if (!vertex)
236     return kFALSE;
237
238   if (vertex->GetNContributors() < 1)
239     return kFALSE;
240
241   ComputeMultiplicity();
242
243   // Variables use to determine the type of trigger :
244   // 0 for minimum bias : CINT1B, CINT1-B or MB1
245   // 1 for muon events : CMUS1B, CMUS1-B or MULow
246   // -1 for everything else
247   TString trigger = fESD->GetFiredTriggerClasses();
248   Int_t triggerClass = -1; 
249
250   if (trigger.Contains("CINT1B") || trigger.Contains("CINT1-B") || trigger.Contains("MB1") || trigger.Contains("CMBAC-B") || trigger.Contains("CMBACS2-B"))
251     triggerClass = 0;
252
253   if (trigger.Contains("CMUS1B") || trigger.Contains("CMUS1-B") || trigger.Contains("MULow"))
254     triggerClass = 1;
255
256   if (triggerClass >= 0)
257     FillHistosESD(triggerClass);
258
259   if (fMCEvent)
260     FillHistosMC();
261
262   return kTRUE;
263 }
264
265
266
267
268 //________________________________________________________________________
269 void AliAnalysisTaskMuonCollisionMultiplicity::FillHistosAOD(Int_t triggerClass)
270 {
271   // Fill histos for AOD events
272   Int_t nTracks = fAOD->GetNumberOfTracks();
273   Int_t nDimuons = fAOD->GetNDimuons();
274   
275   // Fill histos
276   Double_t vertexPosition = fAOD->GetPrimaryVertex()->GetZ();
277   Double_t pileUp = !(fAOD->IsPileupFromSPD());
278   
279   Double_t valuesTrigger[3] = {static_cast<Double_t>(fTrackletMultiplicity), vertexPosition, pileUp};
280   ((THnSparseD *)fTriggerList->At(triggerClass))->Fill(valuesTrigger);
281   
282
283   // Loop on the muons tracks
284   for (Int_t ii = 0; ii < nTracks; ii++){
285     AliAODTrack * aodtrack = dynamic_cast<AliAODTrack*>(fAOD->GetTrack(ii));
286     if(!aodtrack) AliFatal("Not a standard AOD");
287     if (IsUsableMuon(aodtrack))
288       {
289         Double_t matchTrigger = aodtrack->GetMatchTrigger();
290         if (matchTrigger > 1.0)
291           matchTrigger = 1.0;                                            // We don't care what type of trigger it is (low or high pT)
292         
293         Double_t thetaAbs = (180.0 / TMath::Pi()) * TMath::ATan(aodtrack->GetRAtAbsorberEnd()/505.0);
294         Double_t eta = aodtrack->Eta();
295         Double_t p = aodtrack->P();
296         Double_t pT = aodtrack->Pt();
297         // For the p used in the pDCA, we want the mean between the p before the muon go through the absorber (p corrected) and after (p uncorrected)
298         // However p uncorrected is not saved in the AODs
299         // Instead we define p uncorrected as p corrected minus the mean p lost when a muon go through the absorber
300         // p lost for muons (it depends of theta_abs) :
301         // 2.0 < theta_abs < 3.0 ---> 2.98 GeV
302         // 3.0 < theta_abs < 10.0 --->  2.4 GeV
303         // No correction applied otherwise
304         Double_t pDCA = p*aodtrack->DCA();
305         if (2.0 < thetaAbs && thetaAbs < 3.0)
306           pDCA = (p-2.98/2.0) * aodtrack->DCA();
307         if (3.0 < thetaAbs && thetaAbs < 10.0)
308           pDCA = (p-2.4/2.0) * aodtrack->DCA();
309         
310         Double_t valuesMuon[9] = {static_cast<Double_t>(fTrackletMultiplicity), vertexPosition, pileUp, matchTrigger, thetaAbs, eta, pDCA, static_cast<Double_t>(p), static_cast<Double_t>(pT)};
311         ((THnSparseD *)fSingleMuonList->At(triggerClass))->Fill(valuesMuon);
312       }
313   }
314   // Loop on Dimuons
315   for (Int_t ii = 0; ii < nDimuons; ii++)
316     if (fAOD->GetDimuon(ii)->Charge() == 0.0)
317       if (IsUsableMuon(fAOD->GetDimuon(ii)->GetMu(0)))
318         if (IsUsableMuon(fAOD->GetDimuon(ii)->GetMu(1)))
319           {
320             Double_t matchTrigger1 = fAOD->GetDimuon(ii)->GetMu(0)->GetMatchTrigger();
321             if (matchTrigger1 > 0.0)
322               matchTrigger1 = 1.0;
323             Double_t matchTrigger2 = fAOD->GetDimuon(ii)->GetMu(1)->GetMatchTrigger();
324             if (matchTrigger2 > 0.0)
325               matchTrigger2 = 1.0;
326             
327             Double_t nMatchTrigger = matchTrigger1 + matchTrigger2;
328             
329             Double_t thetaAbs1 = (180.0 / TMath::Pi()) * TMath::ATan(fAOD->GetDimuon(ii)->GetMu(0)->GetRAtAbsorberEnd()/505.0);
330             Double_t thetaAbs2 = (180.0 / TMath::Pi()) * TMath::ATan(fAOD->GetDimuon(ii)->GetMu(1)->GetRAtAbsorberEnd()/505.0);
331             Double_t eta1 = fAOD->GetDimuon(ii)->GetMu(0)->Eta();
332             Double_t eta2 = fAOD->GetDimuon(ii)->GetMu(1)->Eta();
333             
334             Double_t p1 = fAOD->GetDimuon(ii)->GetMu(0)->P();
335             Double_t p2 = fAOD->GetDimuon(ii)->GetMu(1)->P();
336             // See the explanation on how the pDCA is computed in the single muon loop
337             Double_t pDCA1 = p1*fAOD->GetDimuon(ii)->GetMu(0)->DCA();
338             if (2.0 < thetaAbs1 && thetaAbs1 < 3.0)
339               pDCA1 = (p1-2.98/2.0) * fAOD->GetDimuon(ii)->GetMu(0)->DCA();
340             if (3.0 < thetaAbs1 && thetaAbs1 < 10.0)
341               pDCA1 = (p1-2.4/2.0) * fAOD->GetDimuon(ii)->GetMu(0)->DCA();
342             Double_t pDCA2 = p2*fAOD->GetDimuon(ii)->GetMu(1)->DCA();
343             if (2.0 < thetaAbs2 && thetaAbs2 < 3.0)
344               pDCA2 = (p2-2.98/2.0) * fAOD->GetDimuon(ii)->GetMu(1)->DCA();
345             if (3.0 < thetaAbs2 && thetaAbs2 < 10.0)
346               pDCA2 = (p2-2.4/2.0) * fAOD->GetDimuon(ii)->GetMu(1)->DCA();
347             
348             Double_t y = fAOD->GetDimuon(ii)->Y();
349             Double_t p = fAOD->GetDimuon(ii)->P();
350             Double_t pT = fAOD->GetDimuon(ii)->Pt();
351             Double_t M = fAOD->GetDimuon(ii)->M();
352             
353             Double_t valuesDimuon[18] = {static_cast<Double_t>(fTrackletMultiplicity), vertexPosition, pileUp, matchTrigger1, matchTrigger2, static_cast<Double_t>(nMatchTrigger), thetaAbs1, thetaAbs2,
354                                          eta1, eta2, pDCA1, pDCA2, static_cast<Double_t>(p1), p2, y, p, pT, M};
355             ((THnSparseD *)fDimuonList->At(triggerClass))->Fill(valuesDimuon);
356           }
357 }
358
359
360
361 void AliAnalysisTaskMuonCollisionMultiplicity::FillHistosESD(Int_t triggerClass)
362 {
363   // Fill the histos for ESD events 
364   Int_t nTracks = fESD->GetNumberOfMuonTracks();
365
366   Double_t vertexPosition = fESD->GetPrimaryVertex()->GetZ();
367   Double_t pileUp = !(fESD->IsPileupFromSPD());
368   
369   Double_t valuesTrigger[3] = {static_cast<Double_t>(fTrackletMultiplicity), vertexPosition, pileUp};
370   ((THnSparseD *)fTriggerList->At(triggerClass))->Fill(valuesTrigger);
371   
372
373   // Loop on the muons tracks
374   for (Int_t ii = 0; ii < nTracks; ii++)
375     if (IsUsableMuon(fESD->GetMuonTrack(ii)))
376       {
377         Double_t matchTrigger1 = fESD->GetMuonTrack(ii)->GetMatchTrigger();
378         if (matchTrigger1 > 1.0)
379           matchTrigger1 = 1.0;                                            // We don't care what type of trigger it is (low or high pT)
380         
381         Double_t thetaAbs1 = (180.0 / TMath::Pi()) * TMath::ATan(fESD->GetMuonTrack(ii)->GetRAtAbsorberEnd()/505.0);
382         Double_t eta1 = fESD->GetMuonTrack(ii)->Eta();
383         Double_t p1 = fESD->GetMuonTrack(ii)->P();
384         Double_t pT1 = fESD->GetMuonTrack(ii)->Pt();
385
386         // For the p used in the pDCA, we want the mean between the p before the muon go through the absorber (p corrected) and after (p uncorrected)
387         // These are respectively AliESDMuonTrack::P() and AliESDMuonTrack::PUncorrected()
388         Double_t pUncor1 = fESD->GetMuonTrack(ii)->PUncorrected();
389         Double_t pDCA1 = (p1+pUncor1) * fESD->GetMuonTrack(ii)->GetDCA() / 2.0;
390         
391         Double_t valuesMuon[9] = {static_cast<Double_t>(fTrackletMultiplicity), vertexPosition, pileUp, matchTrigger1, thetaAbs1, eta1, pDCA1, p1, pT1};
392         ((THnSparseD *)fSingleMuonList->At(triggerClass))->Fill(valuesMuon);
393         
394         // Second loop on muons, to fill the dimuons histos
395         for (Int_t jj = ii+1; jj < nTracks; jj++)
396           if (IsUsableMuon(fESD->GetMuonTrack(jj)))
397             if (fESD->GetMuonTrack(ii)->Charge() + fESD->GetMuonTrack(jj)->Charge() == 0.0)
398               {
399                 Double_t matchTrigger2 = fESD->GetMuonTrack(jj)->GetMatchTrigger();
400                 if (matchTrigger2 > 0.0)
401                   matchTrigger2 = 1.0;
402                 
403                 Double_t nMatchTrigger = matchTrigger1 + matchTrigger2;
404
405                 Double_t thetaAbs2 = (180.0 / TMath::Pi()) * TMath::ATan(fESD->GetMuonTrack(jj)->GetRAtAbsorberEnd()/505.0);
406                 Double_t eta2 = fESD->GetMuonTrack(jj)->Eta();
407                 Double_t p2 = fESD->GetMuonTrack(jj)->P();
408
409                 // For the p used in the pDCA, we want the mean between the p before the muon go through the absorber (p corrected) and after (p uncorrected)
410                 // These are respectively AliESDMuonTrack::P() and AliESDMuonTrack::PUncorrected()
411                 Double_t pUncor2 = fESD->GetMuonTrack(jj)->PUncorrected();
412                 Double_t pDCA2 = (p2+pUncor2) * fESD->GetMuonTrack(jj)->GetDCA() / 2.0;
413                 
414                 // To compute the p, pT and M of the dimuon, we need a TLorentz vector of the dimuon
415                 Double_t E = fESD->GetMuonTrack(ii)->E() + fESD->GetMuonTrack(jj)->E();
416                 Double_t pX = fESD->GetMuonTrack(ii)->Px() + fESD->GetMuonTrack(jj)->Px();
417                 Double_t pY = fESD->GetMuonTrack(ii)->Py() + fESD->GetMuonTrack(jj)->Py();
418                 Double_t pZ = fESD->GetMuonTrack(ii)->Pz() + fESD->GetMuonTrack(jj)->Pz();
419                 TLorentzVector *dimuonVector = new TLorentzVector(pX, pY, pZ, E);
420                 dimuonVector->SetPxPyPzE(pX, pY, pZ, E);
421                 
422                 Double_t y = dimuonVector->Rapidity();
423                 Double_t p = dimuonVector->P();
424                 Double_t pT = TMath::Sqrt(pX*pX + pY*pY);
425                 Double_t M = dimuonVector->M();
426                 
427                 Double_t valuesDimuon[18] = {static_cast<Double_t>(fTrackletMultiplicity), vertexPosition, pileUp, matchTrigger1, matchTrigger2, static_cast<Double_t>(nMatchTrigger), thetaAbs1, thetaAbs2,
428                                              eta1, eta2, pDCA1, pDCA2, p1, p2, y, p, pT, M};
429                 ((THnSparseD *)fDimuonList->At(triggerClass))->Fill(valuesDimuon);
430                 delete dimuonVector;
431               }
432       }
433   
434
435   // Since this is the ESD, we are also going to fill a the V0amp vs multiplicity histos
436   if (triggerClass == 0)
437     {
438       AliESDVZERO *v0 = fESD->GetVZEROData();
439       Float_t multV0 = 0;
440       Float_t multV0A = 0;
441       Float_t multV0C = 0;
442
443       for (Int_t ii = 0; ii < 64; ii++)
444         {
445           multV0 += v0->GetMultiplicity(ii);
446           if (ii < 32)
447             {
448               multV0C += v0->GetMultiplicityV0C(ii);
449               multV0A += v0->GetMultiplicityV0A(ii);
450             }
451         }
452
453       ((TH2D *)fTriggerList->At(2))->Fill(fTrackletMultiplicity, multV0);
454       ((TH2D *)fTriggerList->At(3))->Fill(fTrackletMultiplicity, multV0A);
455       ((TH2D *)fTriggerList->At(4))->Fill(fTrackletMultiplicity, multV0C);
456     }
457
458
459 }
460
461
462 //________________________________________________________________________
463 void AliAnalysisTaskMuonCollisionMultiplicity::FillHistosMC()
464 {
465   // Fill the histo of the correlation between MC tracks and ESD tracks
466
467   Int_t multiplicityGenerated = 0;
468   const AliESDVertex *vertex = fESD->GetPrimaryVertexTracks();
469
470   for (Int_t nn = 0; nn < fMCEvent->GetNumberOfTracks(); nn++)
471     {
472       AliMCParticle *particle = (AliMCParticle *) fMCEvent->GetTrack(nn);
473       Bool_t isGoodMult = kTRUE;
474       
475       if (particle->Particle()->GetStatusCode() != 1)
476         isGoodMult = kFALSE;
477       
478       if (particle->Charge() == 0)
479         isGoodMult = kFALSE;
480
481       if (TMath::Abs(particle->Eta()) > 1.6)
482         isGoodMult = kFALSE;
483
484       // Check if the particle is a pion, kaon, proton, electron or muon
485       if (TMath::Abs(particle->PdgCode()) != 211 && TMath::Abs(particle->PdgCode()) != 321 && TMath::Abs(particle->PdgCode()) != 2212 &&
486           TMath::Abs(particle->PdgCode()) != 11 && TMath::Abs(particle->PdgCode()) != 13)
487         isGoodMult = kFALSE;
488
489       // Check if the distance to vertex is inferior to 1 cm
490       Double_t distanceToVertex = TMath::Sqrt((particle->Xv() - vertex->GetX())*(particle->Xv() - vertex->GetX()) + 
491                                               (particle->Yv() - vertex->GetY())*(particle->Yv() - vertex->GetY()) + 
492                                               (particle->Zv() - vertex->GetZ())*(particle->Zv() - vertex->GetZ()));
493       if (distanceToVertex > 1.0)
494         isGoodMult = kFALSE;
495         
496       if (isGoodMult)
497         multiplicityGenerated += 1;
498     }
499
500   ((TH2D *)fMonteCarloList->At(0))->Fill(multiplicityGenerated, fTrackletMultiplicity);
501 }
502
503
504
505 //________________________________________________________________________
506 void AliAnalysisTaskMuonCollisionMultiplicity::ComputeMultiplicity()
507 {
508   // Compute the collision multiplicity based on AOD or ESD tracklets
509
510   Int_t multiplicity = 0;
511   
512   if (fAOD)
513     {
514       AliAODTracklets *tracklets = fAOD->GetTracklets();
515       Int_t nTracklets = tracklets->GetNumberOfTracklets();
516       for (Int_t nn = 0; nn < nTracklets; nn++)
517         {
518           Double_t theta = tracklets->GetTheta(nn);
519           Double_t eta = -TMath::Log(TMath::Tan(theta/2.0));
520           
521           if (TMath::Abs(eta) < fEtaCut)
522             multiplicity += 1;
523         }
524     }
525
526
527   if (fESD)
528     {
529       // In ESDs, we use the EstimateMultiplicity function
530       Int_t multTracklets = 0;
531       Int_t multTPC = 0;
532       Int_t multITSSA = 0;
533       fESD->EstimateMultiplicity(multTracklets, multTPC, multITSSA, fEtaCut);
534       multiplicity = multTracklets;
535     }
536
537   fTrackletMultiplicity =  multiplicity;
538 }
539
540
541
542 //________________________________________________________________________
543 Bool_t AliAnalysisTaskMuonCollisionMultiplicity::IsUsableMuon(AliAODTrack *track)
544 {
545   // Check if the track is a usable muon track
546   // Cuts applied :
547   // - is it a muon track?
548   // - does it have a pT > 0.0?
549
550   Bool_t isGood = kFALSE;
551
552   if (!track->IsMuonTrack())
553     return isGood;
554
555   if (!(track->Pt() > 0.0))
556     return isGood;
557
558   isGood = kTRUE;
559   return isGood;
560 }
561
562
563
564 //________________________________________________________________________
565 Bool_t AliAnalysisTaskMuonCollisionMultiplicity::IsUsableMuon(AliESDMuonTrack *track)
566 {
567   // Check if the track is a usable muon track
568   // Cuts applied :
569   // - is it a muon track?
570   // - does it have a pT > 0.0?
571
572   Bool_t isGood = kFALSE;
573
574   if (!track->ContainTrackerData())
575     return isGood;
576
577   if (!(track->Pt() > 0.0))
578     return isGood;
579
580   isGood = kTRUE;
581   return isGood;
582 }
583
584
585
586 //________________________________________________________________________
587 void AliAnalysisTaskMuonCollisionMultiplicity::Init()
588 {
589   // Initialize the object
590
591   fTriggerList = new TList();
592   fSingleMuonList = new TList();
593   fDimuonList = new TList();
594   fMonteCarloList = new TList();
595
596   fTriggerList->SetOwner();
597   fSingleMuonList->SetOwner();
598   fDimuonList->SetOwner();
599   fMonteCarloList->SetOwner();
600
601
602
603   // Trigger histos
604   // dimension 0 : multiplicity of the event
605   // dimension 1 : z vertex of the event
606   // dimension 2 : is it an event without pile up (0 for no, 1 for yes)?
607   Int_t nBinsTrigger[3] =       {  150,    60,   2};
608   Double_t minRangeTrigger[3] = {  0.0, -30.0, 0.0};
609   Double_t maxRangeTrigger[3] = {150.0,  30.0, 2.0};
610   THnSparseD *CINT1B = new THnSparseD ("CINT1B", "CINT1B", 3, nBinsTrigger, minRangeTrigger, maxRangeTrigger);
611   THnSparseD *CMUS1B = new THnSparseD ("CMUS1B", "CMUS1B", 3, nBinsTrigger, minRangeTrigger, maxRangeTrigger);
612   CINT1B->Sumw2();
613   CMUS1B->Sumw2();
614
615   TH2D *CompSPDV0 = new TH2D ("CompSPDV0", "CompSPDV0", 150, 0.0, 150.0, 2000, 0.0, 2000.0);
616   CompSPDV0->Sumw2();
617   TH2D *CompSPDV0A = new TH2D ("CompSPDV0A", "CompSPDV0A", 150, 0.0, 150.0, 1000, 0.0, 1000.0);
618   CompSPDV0A->Sumw2();
619   TH2D *CompSPDV0C = new TH2D ("CompSPDV0C", "CompSPDV0C", 150, 0.0, 150.0, 1000, 0.0, 1000.0);
620   CompSPDV0C->Sumw2();
621
622   fTriggerList->AddAt(CINT1B, 0);
623   fTriggerList->AddAt(CMUS1B, 1);
624   fTriggerList->AddAt(CompSPDV0, 2);
625   fTriggerList->AddAt(CompSPDV0A, 3);
626   fTriggerList->AddAt(CompSPDV0C, 4);
627
628
629
630
631   // Muons histos
632   // dimension 0 : multiplicity of the event
633   // dimension 1 : z vertex of the event
634   // dimension 2 : is it an event without pile up (0 for no, 1 for yes)?
635   // dimension 3 : does the muon match the trigger (0 for no, 1 for yes)?
636   // dimension 4 : theta_abs of the muon
637   // dimension 5 : eta of the muon
638   // dimension 6 : p DCA of the muon
639   // dimension 7 : p of the muon
640   // dimension 8 : pT of the muon
641
642   Int_t nBinsMuon[9] =       {  150,    60,   2,   2,  110,   35,   300,   300,  300};
643   Double_t minRangeMuon[9] = {  0.0, -30.0, 0.0, 0.0,  0.0, -5.0,   0.0,   0.0,  0.0};
644   Double_t maxRangeMuon[9] = {150.0,  30.0, 2.0, 2.0, 11.0, -1.5, 300.0, 150.0, 30.0};
645
646   THnSparseD *muonCINT1B = new THnSparseD("muonCINT1B", "muonCINT1B", 9, nBinsMuon, minRangeMuon, maxRangeMuon);
647   THnSparseD *muonCMUS1B = new THnSparseD("muonCMUS1B", "muonCMUS1B", 9, nBinsMuon, minRangeMuon, maxRangeMuon);
648   muonCINT1B->Sumw2();
649   muonCMUS1B->Sumw2();
650
651   fSingleMuonList->AddAt(muonCINT1B, 0);
652   fSingleMuonList->AddAt(muonCMUS1B, 1);
653
654
655   // Dimuons histos
656   // dimension 0  : multiplicity of the event
657   // dimension 1  : z range (0 for no, 1 for yes)?
658   // dimension 2  : is it an event without pile up (0 for no, 1 for yes)?
659   // dimension 3  : does the first muon match the trigger (0 for no, 1 for yes)?
660   // dimension 4  : does the second muon match the trigger (0 for no, 1 for yes)?
661   // dimension 5  : number of muons matching the trigger in the dimuon
662   // dimension 6  : theta_abs of the first muon
663   // dimension 7  : theta_abs of the second muon
664   // dimension 8  : eta of the first muon
665   // dimension 9  : eta of the second muon
666   // dimension 10 : p DCA of the first muon
667   // dimension 11 : p DCA of the second muon
668   // dimension 12 : p of the first muon
669   // dimension 13 : p of the second muon
670   // dimension 14 : y of the dimuon
671   // dimension 15 : p of the dimuon
672   // dimension 16 : pT of the dimuon
673   // dimension 17 : invariant mass of the dimuon
674
675   Int_t nBinsDimuon[18] =       {  150,    60,   2,   2,   2,   3,  110,  110,   35,   35,   200,   200,   300,   300,   35,   300,  300,  375};
676   Double_t minRangeDimuon[18] = {  0.0, -30.0, 0.0, 0.0, 0.0, 0.0,  0.0,  0.0, -5.0, -5.0,   0.0,   0.0,   0.0,   0.0, -5.0,   0.0,  0.0,  0.0};
677   Double_t maxRangeDimuon[18] = {150.0,  30.0, 2.0, 2.0, 2.0, 3.0, 11.0, 11.0, -1.5, -1.5, 600.0, 600.0, 150.0, 150.0, -1.5, 300.0, 30.0, 15.0};
678
679   THnSparseD *dimuonCINT1B = new THnSparseD("dimuonCINT1B", "dimuonCINT1B", 18, nBinsDimuon, minRangeDimuon, maxRangeDimuon);
680   THnSparseD *dimuonCMUS1B = new THnSparseD("dimuonCMUS1B", "dimuonCMUS1B", 18, nBinsDimuon, minRangeDimuon, maxRangeDimuon);
681
682   dimuonCINT1B->Sumw2();
683   dimuonCMUS1B->Sumw2();
684
685   fDimuonList->AddAt(dimuonCINT1B, 0);
686   fDimuonList->AddAt(dimuonCMUS1B, 1);
687
688   // MonteCarlo Histo
689   TH2D *correlGenerReco = new TH2D("correlGenerReco", "correlGenerReco", 250, 0.0, 250.0, 250, 0.0, 250.0);
690   correlGenerReco->GetXaxis()->SetTitle("N ch gener");
691   correlGenerReco->GetYaxis()->SetTitle("N reco tracklets");
692
693   correlGenerReco->Sumw2();
694
695   fMonteCarloList->AddAt(correlGenerReco, 0);
696
697   fIsInit = kTRUE;
698 }
699
700
701
702
703
704 //________________________________________________________________________
705 void AliAnalysisTaskMuonCollisionMultiplicity::Terminate(Option_t */*option*/)
706 {
707 //Terminate analysis
708   
709   fTriggerList = (TList *) GetOutputData(0);
710   fSingleMuonList = (TList *) GetOutputData(1);
711   fDimuonList = (TList *) GetOutputData(2);
712 }