]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/muon/AliAnalysisTaskMuonCollisionMultiplicity.cxx
Coverity fix (an obsolete constructor removed)
[u/mrichter/AliRoot.git] / PWG3 / 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   fZCut(0),
57   fEtaCut(0),
58   fTrackletMultiplicity(0),
59   fTriggerList(0),
60   fSingleMuonList(0),
61   fDimuonList(0),
62   fMonteCarloList(0)
63 {
64   ///Default Constructor
65 }
66
67
68 //___________________________________________________
69 AliAnalysisTaskMuonCollisionMultiplicity::AliAnalysisTaskMuonCollisionMultiplicity(const AliAnalysisTaskMuonCollisionMultiplicity& src)
70   :
71   AliAnalysisTaskSE(),
72   fIsInit(kFALSE),
73   fAOD(0x0),
74   fESD(0x0),
75   fZCut(0),
76   fEtaCut(0),
77   fTrackletMultiplicity(0),
78   fTriggerList(0),
79   fSingleMuonList(0),
80   fDimuonList(0),
81   fMonteCarloList(0)
82 {
83   /// copy ctor
84   src.Copy(*this);
85
86 }
87
88 //___________________________________________________
89 AliAnalysisTaskMuonCollisionMultiplicity& AliAnalysisTaskMuonCollisionMultiplicity::operator=(const AliAnalysisTaskMuonCollisionMultiplicity& src)
90 {
91   /// assignement operator
92   if ( this != &src ) 
93   {
94     src.Copy(*this);
95   }
96   return *this;
97 }
98
99
100
101
102 //___________________________________________________
103 AliAnalysisTaskMuonCollisionMultiplicity::AliAnalysisTaskMuonCollisionMultiplicity(const Char_t *name)
104   :
105   AliAnalysisTaskSE(name),
106   fIsInit(kFALSE),
107   fAOD(0x0),
108   fESD(0x0),
109   fZCut(0),
110   fEtaCut(0),
111   fTrackletMultiplicity(0),
112   fTriggerList(0),
113   fSingleMuonList(0),
114   fDimuonList(0),
115   fMonteCarloList(0)
116 {
117   // Define Inputs and outputs
118   //DefineInput(0, TChain::Class());
119   //DefineInput(1, TChain::Class());
120   DefineOutput(0, TList::Class());
121   DefineOutput(1, TList::Class());
122   DefineOutput(2, TList::Class());
123   DefineOutput(3, TList::Class());
124   DefineOutput(4, TList::Class());
125 }
126
127
128 //______________________________________________________________________________
129 AliAnalysisTaskMuonCollisionMultiplicity::~AliAnalysisTaskMuonCollisionMultiplicity()
130 {
131 // Destructor.
132   delete fAOD;
133   delete fESD;
134   delete fTriggerList;
135   delete fSingleMuonList;
136   delete fDimuonList;
137   delete fMonteCarloList;
138 }
139
140
141
142 //________________________________________________________________________
143 void AliAnalysisTaskMuonCollisionMultiplicity::UserCreateOutputObjects()
144 {
145   if (!fIsInit)
146     Init();
147   OpenFile(0);
148 }
149
150
151
152
153 //________________________________________________________________________
154 void AliAnalysisTaskMuonCollisionMultiplicity::UserExec(Option_t */*option*/)
155 {
156   fAOD = 0x0;
157   fESD = 0x0;
158
159   if (!fIsInit)
160     Init();
161
162   fAOD = dynamic_cast<AliAODEvent *> (InputEvent());
163   if (!fAOD)
164     fESD = dynamic_cast<AliESDEvent *> (InputEvent());
165   
166  
167   if (fAOD)
168     CheckEventAOD();
169
170   if (fESD)
171     CheckEventESD();
172
173   PostData(1, fTriggerList);
174   PostData(2, fSingleMuonList);
175   PostData(3, fDimuonList);
176   PostData(4, fMonteCarloList);
177 }
178
179
180
181 //________________________________________________________________________
182 void AliAnalysisTaskMuonCollisionMultiplicity::NotifyRun()
183 {
184   
185 }
186
187
188 //________________________________________________________________________
189 void AliAnalysisTaskMuonCollisionMultiplicity::FinishTaskOutput()
190 {
191   
192 }
193
194
195 //__________________________________________________________________________
196 Bool_t AliAnalysisTaskMuonCollisionMultiplicity::CheckEventAOD()
197 {
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   if (fTrackletMultiplicity < 1)
208     return kFALSE;
209   
210   // Variables use to determine the type of trigger :
211   // 0 for minimum bias : CINT1B, CINT1-B or MB1
212   // 1 for muon events : CMUS1B, CMUS1-B or MULow
213   // -1 for everything else
214   TString trigger = fAOD->GetFiredTriggerClasses();
215   Int_t triggerClass = -1; 
216
217   if (trigger.Contains("CINT1B") || trigger.Contains("CINT1-B") || trigger.Contains("MB1"))
218     triggerClass = 0;
219
220   if (trigger.Contains("CMUS1B") || trigger.Contains("CMUS1-B") || trigger.Contains("MULow"))
221     triggerClass = 1;
222
223   if (triggerClass >= 0)
224     FillHistosAOD(triggerClass);
225   
226   return kTRUE;
227 }
228
229
230
231 //__________________________________________________________________________
232 Bool_t AliAnalysisTaskMuonCollisionMultiplicity::CheckEventESD()
233 {
234   const AliESDVertex *vertex = fESD->GetPrimaryVertex();
235   
236   if (!vertex)
237     return kFALSE;
238
239   if (vertex->GetNContributors() < 1)
240     return kFALSE;
241
242   ComputeMultiplicity();
243   if (fTrackletMultiplicity < 1)
244     return kFALSE;
245   
246   // Variables use to determine the type of trigger :
247   // 0 for minimum bias : CINT1B, CINT1-B or MB1
248   // 1 for muon events : CMUS1B, CMUS1-B or MULow
249   // -1 for everything else
250   TString trigger = fESD->GetFiredTriggerClasses();
251   Int_t triggerClass = -1; 
252
253   if (trigger.Contains("CINT1B") || trigger.Contains("CINT1-B") || trigger.Contains("MB1") || trigger.Contains("CMBAC-B") || trigger.Contains("CMBACS2-B"))
254     triggerClass = 0;
255
256   if (trigger.Contains("CMUS1B") || trigger.Contains("CMUS1-B") || trigger.Contains("MULow"))
257     triggerClass = 1;
258
259   if (triggerClass >= 0)
260     FillHistosESD(triggerClass);
261
262   if (fMCEvent)
263     FillHistosMC();
264
265   return kTRUE;
266 }
267
268
269
270
271 //________________________________________________________________________
272 void AliAnalysisTaskMuonCollisionMultiplicity::FillHistosAOD(Int_t triggerClass)
273 {
274   // Fill histos for AOD events
275   Int_t nTracks = fAOD->GetNTracks();
276   Int_t nDimuons = fAOD->GetNDimuons();
277   
278   // Fill histos
279   Double_t vertexCut = (TMath::Abs(fAOD->GetPrimaryVertex()->GetZ()) < fZCut);
280   Double_t pileUp = !(fAOD->IsPileupFromSPD());
281   
282   Double_t valuesTrigger[3] = {fTrackletMultiplicity, vertexCut, pileUp};
283   ((THnSparseD *)fTriggerList->At(triggerClass))->Fill(valuesTrigger);
284   
285   // Loop on the muons tracks
286   for (Int_t ii = 0; ii < nTracks; ii++)
287     if (IsUsableMuon(fAOD->GetTrack(ii)))
288       {
289         Double_t matchTrigger = fAOD->GetTrack(ii)->GetMatchTrigger();
290         if (matchTrigger > 1.0)
291           matchTrigger = 1.0;                                            // We don't care what type of trigger it is
292         
293         Double_t thetaAbs = (180.0 / TMath::Pi()) * TMath::ATan(fAOD->GetTrack(ii)->GetRAtAbsorberEnd()/505.0);
294         Double_t eta = fAOD->GetTrack(ii)->Eta();
295         Double_t dcaX = fAOD->GetTrack(ii)->XAtDCA();
296         Double_t dcaY = fAOD->GetTrack(ii)->YAtDCA();
297         Double_t p = fAOD->GetTrack(ii)->P();
298         Double_t pT = fAOD->GetTrack(ii)->Pt();
299         Double_t pDCA = p*TMath::Sqrt(dcaX*dcaX + dcaY*dcaY);
300         
301         Double_t valuesMuon[11] = {fTrackletMultiplicity, vertexCut, pileUp, matchTrigger, thetaAbs, eta, p*dcaX, p*dcaY, pDCA, p, pT};
302         ((THnSparseD *)fSingleMuonList->At(triggerClass))->Fill(valuesMuon);
303       }
304   
305   // Loop on Dimuons
306   for (Int_t ii = 0; ii < nDimuons; ii++)
307     if (fAOD->GetDimuon(ii)->Charge() == 0.0)
308       if (IsUsableMuon(fAOD->GetDimuon(ii)->GetMu(0)))
309         if (IsUsableMuon(fAOD->GetDimuon(ii)->GetMu(1)))
310           {
311             Double_t matchTrigger1 = fAOD->GetDimuon(ii)->GetMu(0)->GetMatchTrigger();
312             if (matchTrigger1 > 0.0)
313               matchTrigger1 = 1.0;
314             Double_t matchTrigger2 = fAOD->GetDimuon(ii)->GetMu(1)->GetMatchTrigger();
315             if (matchTrigger2 > 0.0)
316               matchTrigger2 = 1.0;
317             
318             Double_t nMatchTrigger = matchTrigger1 + matchTrigger2;
319             
320             Double_t thetaAbs1 = (180.0 / TMath::Pi()) * TMath::ATan(fAOD->GetDimuon(ii)->GetMu(0)->GetRAtAbsorberEnd()/505.0);
321             Double_t thetaAbs2 = (180.0 / TMath::Pi()) * TMath::ATan(fAOD->GetDimuon(ii)->GetMu(1)->GetRAtAbsorberEnd()/505.0);
322             Double_t eta1 = fAOD->GetDimuon(ii)->GetMu(0)->Eta();
323             Double_t eta2 = fAOD->GetDimuon(ii)->GetMu(1)->Eta();
324             Double_t dcaX1 = fAOD->GetDimuon(ii)->GetMu(0)->XAtDCA();
325             Double_t dcaY1 = fAOD->GetDimuon(ii)->GetMu(0)->YAtDCA();
326             Double_t dcaX2 = fAOD->GetDimuon(ii)->GetMu(1)->XAtDCA();
327             Double_t dcaY2 = fAOD->GetDimuon(ii)->GetMu(1)->YAtDCA();
328             
329             Double_t p1 = fAOD->GetDimuon(ii)->GetMu(0)->P();
330             Double_t p2 = fAOD->GetDimuon(ii)->GetMu(1)->P();
331             
332             Double_t p = fAOD->GetDimuon(ii)->P();
333             Double_t pT = fAOD->GetDimuon(ii)->Pt();
334             Double_t M = fAOD->GetDimuon(ii)->M();
335             
336             Double_t valuesDimuon[19] = {fTrackletMultiplicity, vertexCut, pileUp, matchTrigger1, matchTrigger2, nMatchTrigger, thetaAbs1, thetaAbs2,
337                                          eta1, eta2, p1*dcaX1, p1*dcaY1, p2*dcaX2, p2*dcaY2, p1, p2, p, pT, M};
338             ((THnSparseD *)fDimuonList->At(triggerClass))->Fill(valuesDimuon);
339           }
340 }
341
342
343
344 void AliAnalysisTaskMuonCollisionMultiplicity::FillHistosESD(Int_t triggerClass)
345 {
346   // Fill the histos for ESD events 
347   // This is not yet tested (in particular the Dimuon part)
348
349   Int_t nTracks = fESD->GetNumberOfMuonTracks();
350
351   Double_t vertexCut = (TMath::Abs(fESD->GetPrimaryVertex()->GetZ()) < fZCut);
352   Double_t pileUp = !(fESD->IsPileupFromSPD());
353   
354   Double_t valuesTrigger[3] = {fTrackletMultiplicity, vertexCut, pileUp};
355   ((THnSparseD *)fTriggerList->At(triggerClass))->Fill(valuesTrigger);
356   
357   // Loop on the muons tracks
358   for (Int_t ii = 0; ii < nTracks; ii++)
359     if (IsUsableMuon(fESD->GetMuonTrack(ii)))
360       {
361         Double_t matchTrigger1 = fESD->GetMuonTrack(ii)->GetMatchTrigger();
362         if (matchTrigger1 > 1.0)
363           matchTrigger1 = 1.0;                                            // We don't care what type of trigger it is
364         
365         Double_t thetaAbs1 = (180.0 / TMath::Pi()) * TMath::ATan(fESD->GetMuonTrack(ii)->GetRAtAbsorberEnd()/505.0);
366         Double_t eta1 = fESD->GetMuonTrack(ii)->Eta();
367         Double_t dcaX1 = fESD->GetMuonTrack(ii)->GetNonBendingCoorAtDCA();
368         Double_t dcaY1 = fESD->GetMuonTrack(ii)->GetBendingCoorAtDCA();
369         Double_t p1 = fESD->GetMuonTrack(ii)->P();
370         Double_t pUncor1 = fESD->GetMuonTrack(ii)->PUncorrected();
371         Double_t pT1 = fESD->GetMuonTrack(ii)->Pt();
372
373         // The pDCA is computed with the average of p and pUncor
374         Double_t pDCAX1 = (p1+pUncor1) * dcaX1 / 2.0;
375         Double_t pDCAY1 = (p1+pUncor1) * dcaY1 / 2.0;
376         Double_t pDCA1 = (p1+pUncor1) * TMath::Sqrt(dcaX1*dcaX1 + dcaY1*dcaY1) / 2.0;
377         
378         Double_t valuesMuon[11] = {fTrackletMultiplicity, vertexCut, pileUp, matchTrigger1, thetaAbs1, eta1, pDCAX1, pDCAY1, pDCA1, p1, pT1};
379         ((THnSparseD *)fSingleMuonList->At(triggerClass))->Fill(valuesMuon);
380         
381         // Second loop on muons, to fill the dimuons histos
382         for (Int_t jj = ii+1; jj < nTracks; jj++)
383           if (IsUsableMuon(fESD->GetMuonTrack(jj)))
384             if (fESD->GetMuonTrack(ii)->Charge() + fESD->GetMuonTrack(jj)->Charge() == 0.0)
385               {
386                 Double_t matchTrigger2 = fESD->GetMuonTrack(jj)->GetMatchTrigger();
387                 if (matchTrigger2 > 0.0)
388                   matchTrigger2 = 1.0;
389                 
390                 Double_t nMatchTrigger = matchTrigger1 + matchTrigger2;
391
392                 Double_t thetaAbs2 = (180.0 / TMath::Pi()) * TMath::ATan(fESD->GetMuonTrack(jj)->GetRAtAbsorberEnd()/505.0);
393                 Double_t eta2 = fESD->GetMuonTrack(jj)->Eta();
394                 Double_t dcaX2 = fESD->GetMuonTrack(jj)->GetNonBendingCoorAtDCA();
395                 Double_t dcaY2 = fESD->GetMuonTrack(jj)->GetBendingCoorAtDCA();
396                 Double_t p2 = fESD->GetMuonTrack(jj)->P();
397                 Double_t pUncor2 = fESD->GetMuonTrack(jj)->PUncorrected();
398
399                 // The pDCA is computed with the average of p and pUncor
400                 Double_t pDCAX2 = (p2+pUncor2) * dcaX2/2.0;
401                 Double_t pDCAY2 = (p2+pUncor2) * dcaY2/2.0;
402                 
403                 // To compute the p, pT and M of the dimuon, we need a TLorentz vector of the dimuon
404                 Double_t E = fESD->GetMuonTrack(ii)->E() + fESD->GetMuonTrack(jj)->E();
405                 Double_t pX = fESD->GetMuonTrack(ii)->Px() + fESD->GetMuonTrack(jj)->Px();
406                 Double_t pY = fESD->GetMuonTrack(ii)->Py() + fESD->GetMuonTrack(jj)->Py();
407                 Double_t pZ = fESD->GetMuonTrack(ii)->Pz() + fESD->GetMuonTrack(jj)->Pz();
408                 TLorentzVector *dimuonVector = new TLorentzVector(pX, pY, pZ, E);
409                 dimuonVector->SetPxPyPzE(pX, pY, pZ, E);
410                 
411                 Double_t p = dimuonVector->P();
412                 Double_t pT = TMath::Sqrt(pX*pX + pY*pY);
413                 Double_t M = dimuonVector->M();
414                 
415                 Double_t valuesDimuon[19] = {fTrackletMultiplicity, vertexCut, pileUp, matchTrigger1, matchTrigger2, nMatchTrigger, thetaAbs1, thetaAbs2,
416                                              eta1, eta2, pDCAX1, pDCAY1, pDCAX2, pDCAY2, p1, p2, p, pT, M};
417                 ((THnSparseD *)fDimuonList->At(triggerClass))->Fill(valuesDimuon);
418                 delete dimuonVector;
419               }
420       }
421   
422 }
423
424
425 //________________________________________________________________________
426 void AliAnalysisTaskMuonCollisionMultiplicity::FillHistosMC()
427 {
428   // Fill the histo of the correlation between MC tracks and ESD tracks
429
430   Int_t multiplicityGenerated = 0;
431   const AliESDVertex *vertex = fESD->GetPrimaryVertexTracks();
432
433   for (Int_t nn = 0; nn < fMCEvent->GetNumberOfTracks(); nn++)
434     {
435       AliMCParticle *particle = (AliMCParticle *) fMCEvent->GetTrack(nn);
436       Bool_t isGoodMult = kTRUE;
437       
438       if (particle->Particle()->GetStatusCode() != 1)
439         isGoodMult = kFALSE;
440       
441       if (particle->Charge() == 0)
442         isGoodMult = kFALSE;
443
444       if (TMath::Abs(particle->Eta()) > 1.0)
445         isGoodMult = kFALSE;
446
447       // Check if the particle is a pion, kaon, proton, electron or muon
448       if (TMath::Abs(particle->PdgCode()) != 211 && TMath::Abs(particle->PdgCode()) != 321 && TMath::Abs(particle->PdgCode()) != 2212 &&
449           TMath::Abs(particle->PdgCode()) != 11 && TMath::Abs(particle->PdgCode()) != 13)
450         isGoodMult = kFALSE;
451
452       // Check if the distance to vertex is inferior to 1 cm
453       Double_t distanceToVertex = TMath::Sqrt((particle->Xv() - vertex->GetXv())*(particle->Xv() - vertex->GetXv()) + 
454                                               (particle->Yv() - vertex->GetYv())*(particle->Yv() - vertex->GetYv()) + 
455                                               (particle->Zv() - vertex->GetZv())*(particle->Zv() - vertex->GetZv()));
456       if (distanceToVertex > 1.0)
457         isGoodMult = kFALSE;
458         
459       if (isGoodMult)
460         multiplicityGenerated += 1;
461     }
462
463   ((TH2D *)fMonteCarloList->At(0))->Fill(multiplicityGenerated, fTrackletMultiplicity);
464 }
465
466
467
468 //________________________________________________________________________
469 void AliAnalysisTaskMuonCollisionMultiplicity::ComputeMultiplicity()
470 {
471   // Compute the collision multiplicity based on AOD or ESD tracklets
472
473   Int_t multiplicity = 0;
474   
475   if (fAOD)
476     {
477       AliAODTracklets *tracklets = fAOD->GetTracklets();
478       Int_t nTracklets = tracklets->GetNumberOfTracklets();
479       for (Int_t nn = 0; nn < nTracklets; nn++)
480         {
481           Double_t theta = tracklets->GetTheta(nn);
482           Double_t eta = -TMath::Log(TMath::Tan(theta/2.0));
483           
484           if (TMath::Abs(eta) < fEtaCut)
485             multiplicity += 1;
486         }
487     }
488
489
490   else
491     for (Int_t nn = 0; nn < fESD->GetMultiplicity()->GetNumberOfTracklets(); nn++)
492       if (TMath::Abs(fESD->GetMultiplicity()->GetEta(nn)) < fEtaCut)
493         multiplicity += 1;
494       
495     
496
497   fTrackletMultiplicity =  multiplicity;
498 }
499
500
501
502 //________________________________________________________________________
503 Bool_t AliAnalysisTaskMuonCollisionMultiplicity::IsUsableMuon(AliAODTrack *track)
504 {
505   // Check if the track is a usable muon track
506   // Cuts applied :
507   // - is it a muon track?
508   // - does it have a pT > 0.0?
509
510   Bool_t isGood = kFALSE;
511
512   if (!track->IsMuonTrack())
513     return isGood;
514
515   if (!(track->Pt() > 0.0))
516     return isGood;
517
518   isGood = kTRUE;
519   return isGood;
520 }
521
522
523
524 //________________________________________________________________________
525 Bool_t AliAnalysisTaskMuonCollisionMultiplicity::IsUsableMuon(AliESDMuonTrack *track)
526 {
527   // Check if the track is a usable muon track
528   // Cuts applied :
529   // - is it a muon track?
530   // - does it have a pT > 0.0?
531
532   Bool_t isGood = kFALSE;
533
534   if (!track->ContainTrackerData())
535     return isGood;
536
537   if (!(track->Pt() > 0.0))
538     return isGood;
539
540   isGood = kTRUE;
541   return isGood;
542 }
543
544
545
546 //________________________________________________________________________
547 void AliAnalysisTaskMuonCollisionMultiplicity::Init()
548 {
549   // Initialize the object
550
551   fTriggerList = new TList();
552   fSingleMuonList = new TList();
553   fDimuonList = new TList();
554   fMonteCarloList = new TList();
555
556   fTriggerList->SetOwner();
557   fSingleMuonList->SetOwner();
558   fDimuonList->SetOwner();
559   fMonteCarloList->SetOwner();
560
561
562
563
564   // Trigger histos
565   // dimension 0 : multiplicity of the event
566   // dimension 1 : is the vertex in the z range (0 for no, 1 for yes)?
567   // dimension 2 : is it an event without pile up (0 for no, 1 for yes)?
568   Int_t nBinsTrigger[3] =       {  150,   2,   2};
569   Double_t minRangeTrigger[3] = {  0.0, 0.0, 0.0};
570   Double_t maxRangeTrigger[3] = {150.0, 2.0, 2.0};
571   THnSparseD *CINT1B = new THnSparseD ("CINT1B", "CINT1B", 3, nBinsTrigger, minRangeTrigger, maxRangeTrigger);
572   THnSparseD *CMUS1B = new THnSparseD ("CMUS1B", "CMUS1B", 3, nBinsTrigger, minRangeTrigger, maxRangeTrigger);
573   CINT1B->Sumw2();
574   CMUS1B->Sumw2();
575
576   fTriggerList->AddAt(CINT1B, 0);
577   fTriggerList->AddAt(CMUS1B, 1);
578
579
580
581
582
583   // Muons histos
584   // dimension 0  : multiplicity of the event
585   // dimension 1  : is the vertex in the z range (0 for no, 1 for yes)?
586   // dimension 2  : is it an event without pile up (0 for no, 1 for yes)?
587   // dimension 3  : does the muon match the trigger (0 for no, 1 for yes)?
588   // dimension 4  : theta_abs of the muon
589   // dimension 5  : eta of the muon
590   // dimension 6  : p DCA_x of the muon
591   // dimension 7  : p DCA_y of the muon
592   // dimension 8  : p DCA of the muon
593   // dimension 9  : p of the muon
594   // dimension 10 : pT of the muon
595
596   Int_t nBinsMuon[11] =       {  150,   2,   2,   2,  110,   35,    150,    150,   150,   500,  300};
597   Double_t minRangeMuon[11] = {  0.0, 0.0, 0.0, 0.0,  0.0, -5.0, -300.0, -300.0,   0.0,   0.0,  0.0};
598   Double_t maxRangeMuon[11] = {150.0, 2.0, 2.0, 2.0, 11.0, -1.5,  300.0,  300.0, 450.0, 100.0, 30.0};
599
600   THnSparseD *muonCINT1B = new THnSparseD("muonCINT1B", "muonCINT1B", 11, nBinsMuon, minRangeMuon, maxRangeMuon);
601   THnSparseD *muonCMUS1B = new THnSparseD("muonCMUS1B", "muonCMUS1B", 11, nBinsMuon, minRangeMuon, maxRangeMuon);
602   muonCINT1B->Sumw2();
603   muonCMUS1B->Sumw2();
604
605   fSingleMuonList->AddAt(muonCINT1B, 0);
606   fSingleMuonList->AddAt(muonCMUS1B, 1);
607
608
609   // Dimuons histos
610   // dimension 0  : multiplicity of the event
611   // dimension 1  : is the vertex in the z range (0 for no, 1 for yes)?
612   // dimension 2  : is it an event without pile up (0 for no, 1 for yes)?
613   // dimension 3  : does the first muon match the trigger (0 for no, 1 for yes)?
614   // dimension 4  : does the second muon match the trigger (0 for no, 1 for yes)?
615   // dimension 5  : number of muons matching the trigger in the dimuon
616   // dimension 6  : theta_abs of the first muon
617   // dimension 7  : theta_abs of the second muon
618   // dimension 8  : eta of the first muon
619   // dimension 9  : eta of the second muon
620   // dimension 10 : p DCA_x of the first muon
621   // dimension 11 : p DCA_y of the first muon
622   // dimension 12 : p DCA_x of the second muon
623   // dimension 13 : p DCA_y of the second muon
624   // dimension 14 : p of the first muon
625   // dimension 15 : p of the second muon
626   // dimension 16 : p of the dimuon
627   // dimension 17 : pT of the dimuon
628   // dimension 18 : invariant mass of the dimuon
629
630   Int_t nBinsDimuon[19] =       {  150,   2,   2,   2,   2,   3,  110,  110,   35,   35,    150,    150,    150,    150,   500,   500,   500,  300,  375};
631   Double_t minRangeDimuon[19] = {  0.0, 0.0, 0.0, 0.0, 0.0, 0.0,  0.0,  0.0, -5.0, -5.0, -300.0, -300.0, -300.0, -300.0,   0.0,   0.0,   0.0,  0.0,  0.0};
632   Double_t maxRangeDimuon[19] = {150.0, 2.0, 2.0, 2.0, 2.0, 3.0, 11.0, 11.0, -1.5, -1.5,  300.0,  300.0,  300.0,  300.0, 100.0, 100.0, 100.0, 30.0, 15.0};
633
634   THnSparseD *dimuonCINT1B = new THnSparseD("dimuonCINT1B", "dimuonCINT1B", 19, nBinsDimuon, minRangeDimuon, maxRangeDimuon);
635   THnSparseD *dimuonCMUS1B = new THnSparseD("dimuonCMUS1B", "dimuonCMUS1B", 19, nBinsDimuon, minRangeDimuon, maxRangeDimuon);
636
637   dimuonCINT1B->Sumw2();
638   dimuonCMUS1B->Sumw2();
639
640   fDimuonList->AddAt(dimuonCINT1B, 0);
641   fDimuonList->AddAt(dimuonCMUS1B, 1);
642
643   // MonteCarlo Histo
644   TH2D *correlGenerReco = new TH2D("correlGenerReco", "correlGenerReco", 250, 0.0, 250.0, 250, 0.0, 250.0);
645   correlGenerReco->GetXaxis()->SetTitle("N ch gener");
646   correlGenerReco->GetYaxis()->SetTitle("N reco tracklets");
647
648   correlGenerReco->Sumw2();
649
650   fMonteCarloList->AddAt(correlGenerReco, 0);
651
652   fIsInit = kTRUE;
653 }
654
655
656
657
658
659 //________________________________________________________________________
660 void AliAnalysisTaskMuonCollisionMultiplicity::Terminate(Option_t */*option*/)
661 {
662 //Terminate analysis
663   
664   fTriggerList = (TList *) GetOutputData(0);
665   fSingleMuonList = (TList *) GetOutputData(1);
666   fDimuonList = (TList *) GetOutputData(2);
667 }