]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/muon/AliAnalysisMuMuSingle.cxx
Merge branch 'TPCdev' of https://git.cern.ch/reps/AliRoot into TPCdev
[u/mrichter/AliRoot.git] / PWG / muon / AliAnalysisMuMuSingle.cxx
1 #include "AliAnalysisMuMuSingle.h"
2
3 /**
4  *
5  * \ingroup pwg-muon-mumu
6  *
7  * \class AliAnalysisMuMuSingle
8  *
9  * Histogramming of single muon tracks. Mostly to get control plots for
10  * the AliAnalysisMuMuMinv sub-analysis, with respect to track cuts used,
11  * like Rabs, p x DCA, etc...
12  *
13  */
14
15 #include "TH2F.h"
16 #include "AliMuonTrackCuts.h"
17 #include "AliAnalysisMuonUtility.h"
18 #include "TMath.h"
19 #include "AliLog.h"
20 #include "AliVParticle.h"
21 #include "TLorentzVector.h"
22 #include "AliAnalysisMuMuCutCombination.h"
23 #include "AliAnalysisMuMuCutRegistry.h"
24 #include "AliMergeableCollection.h"
25
26 ClassImp(AliAnalysisMuMuSingle)
27
28 //_____________________________________________________________________________
29 AliAnalysisMuMuSingle::AliAnalysisMuMuSingle()
30 : AliAnalysisMuMuBase(),
31 fMuonTrackCuts(0x0),
32 fShouldSeparatePlusAndMinus(kFALSE),
33 fAccEffHisto(0x0)
34 {
35   /// ctor
36 }
37
38 //_____________________________________________________________________________
39 AliAnalysisMuMuSingle::~AliAnalysisMuMuSingle()
40 {
41   /// dtor
42   delete fMuonTrackCuts;
43   delete fAccEffHisto;
44 }
45
46
47 //_____________________________________________________________________________
48 void
49 AliAnalysisMuMuSingle::CreateTrackHisto(const char* eventSelection,
50                                         const char* triggerClassName,
51                                         const char* centrality,
52                                         const char* hname, const char* htitle,
53                                         Int_t nbinsx, Double_t xmin, Double_t xmax,
54                                         Int_t nbinsy, Double_t ymin, Double_t ymax,
55                                         Bool_t separatePlusAndMinus) const
56 {
57   /// Append histograms for single track to our histogram collection
58   
59   if ( IsHistogramDisabled(hname) ) return;
60   
61   if ( separatePlusAndMinus )
62   {
63     const char* suffix[] = { "Plus", "Minus" };
64     const char* symbol[] = { "+", "-" };
65     
66     for ( Int_t i = 0; i < 2; ++i )
67     {
68       TString shtitle(htitle);
69       TString shname(hname);
70       
71       shtitle.ReplaceAll("#mu",Form("#mu^{%s}",symbol[i]));
72       
73       shname += suffix[i];
74       
75       CreateTrackHistos(kHistoForData | kHistoForMCInput,eventSelection,triggerClassName,centrality,shname.Data(),shtitle.Data(),
76                         nbinsx,xmin,xmax,nbinsy,ymin,ymax);
77     }
78   }
79   else
80   {
81     CreateTrackHistos(kHistoForData | kHistoForMCInput,eventSelection,triggerClassName,centrality,hname,htitle,
82                 nbinsx,xmin,xmax,nbinsy,ymin,ymax);
83   }
84 }
85
86 //_____________________________________________________________________________
87 Bool_t AliAnalysisMuMuSingle::IsPDCAOK(const AliVParticle& part)
88 {
89   UInt_t selectionMask = MuonTrackCuts() ? MuonTrackCuts()->GetSelectionMask(&part) : 0;
90
91   return ( ( selectionMask & AliMuonTrackCuts::kMuPdca ) == AliMuonTrackCuts::kMuPdca );
92 }
93
94 //_____________________________________________________________________________
95 Bool_t AliAnalysisMuMuSingle::IsRabsOK(const AliVParticle& part) const
96 {
97   Double_t thetaAbsEndDeg = AliAnalysisMuonUtility::GetThetaAbsDeg(&part);
98   
99   return ( thetaAbsEndDeg > 2. && thetaAbsEndDeg < 10. );
100 }
101
102 //_____________________________________________________________________________
103 Bool_t AliAnalysisMuMuSingle::IsEtaInRange(const AliVParticle& part) const
104 {
105   return (part.Eta() > -4.0 && part.Eta() < -2.5);
106 }
107
108 //_____________________________________________________________________________
109 void AliAnalysisMuMuSingle::SetRun(const AliInputEventHandler* eventHandler)
110 {
111   MuonTrackCuts()->SetRun(eventHandler);
112 }
113
114 //_____________________________________________________________________________
115 Int_t AliAnalysisMuMuSingle::EAGetNumberOfMuonTracks() const
116 {
117   // Get the number of muon tracks *that are not ghosts*
118   
119   Int_t ntracks = AliAnalysisMuonUtility::GetNTracks(Event());
120   
121   for ( Int_t i = 0; i < ntracks; ++i )
122   {
123     AliVParticle* track = AliAnalysisMuonUtility::GetTrack(i,Event());
124     if (AliAnalysisMuonUtility::IsMuonGhost(track)) --ntracks;
125   }
126   
127   return ntracks;
128 }
129
130 ////_____________________________________________________________________________
131 //Int_t AliAnalysisMuMuSingle::EAGetNumberOfSelectMuonTracks() const
132 //{
133 //  // Get the number of "very good" muon tracks :
134 //  // Rabs + DCA + pT > 1.5 Gev/C
135 //  
136 //  Int_t nTracks = AliAnalysisMuonUtility::GetNTracks(Event());
137 //  
138 //  UInt_t check = kAll | kMatched | kRabs | kDCA | kEta | kPt1dot5;
139 //  
140 //  Int_t nGood(0);
141 //  
142 //  for ( Int_t i = 0; i < nTracks; ++i )
143 //  {
144 //    ULong64_t m = GetTrackMask(i);
145 //    if ( ( m & check ) == check )
146 //    {
147 //      ++nGood;
148 //    }
149 //  }
150 //  return nGood;
151 //}
152
153 //_____________________________________________________________________________
154 Double_t AliAnalysisMuMuSingle::EAGetTrackDCA(const AliVParticle& track) const
155 {
156   // Get track DCA
157   
158   Double_t xdca = AliAnalysisMuonUtility::GetXatDCA(&track);
159   Double_t ydca = AliAnalysisMuonUtility::GetYatDCA(&track);
160   
161   return TMath::Sqrt(xdca*xdca+ydca*ydca);
162 }
163
164 //_____________________________________________________________________________
165 void AliAnalysisMuMuSingle::DefineHistogramCollection(const char* eventSelection,
166                                                       const char* triggerClassName,
167                                                       const char* centrality)
168 {
169   /// Actually create the histograms for phyics/triggerClassName
170  
171   if ( Histo(eventSelection,triggerClassName,centrality,"AliAnalysisMuMuSingle") )
172   {
173     return;
174   }
175
176   AliAnalysisMuMuBase::EDataType dt = AliAnalysisMuMuBase::kHistoForData;
177   
178   // dummy histogram to signal that we already defined all our histograms (see above)
179   CreateEventHistos(dt,eventSelection,triggerClassName,centrality,"AliAnalysisMuMuSingle","Dummy semaphore",1,0,1);
180   
181   Double_t ptMin = 0;
182   Double_t ptMax = 12*3;
183   Int_t nbinsPt = GetNbins(ptMin,ptMax,0.5);
184   Double_t pMin = 0;
185   Double_t pMax = 100*3;
186   Int_t nbinsP = GetNbins(pMin,pMax,2.0);
187   Double_t etaMin = -5;
188   Double_t etaMax = -2;
189   Int_t nbinsEta = GetNbins(etaMin,etaMax,0.05);
190   
191   Double_t rapidityMin = -5;
192   Double_t rapidityMax = -2;
193   Int_t nbinsRapidity = GetNbins(rapidityMin,rapidityMax,0.05);
194   
195   Double_t phiMin = -TMath::Pi();
196   Double_t phiMax = TMath::Pi();
197   Int_t nbinsPhi = GetNbins(phiMin,phiMax,0.05);
198   
199   CreateTrackHisto(eventSelection,triggerClassName,centrality,"Chi2MatchTrigger","Chi2 Match Trigger",72,0,72);
200   
201   CreateTrackHisto(eventSelection,triggerClassName,centrality,"EtaRapidityMu", "Eta distribution vs Rapidity for #mu", nbinsRapidity,rapidityMin,rapidityMax,nbinsEta,etaMin,etaMax, fShouldSeparatePlusAndMinus);
202   
203   CreateTrackHisto(eventSelection,triggerClassName,centrality,"PtEtaMu", "P_{T} distribution vs Eta for #mu", nbinsEta,etaMin,etaMax, nbinsPt,ptMin,ptMax,fShouldSeparatePlusAndMinus);
204   
205   CreateTrackHisto(eventSelection,triggerClassName,centrality,"PtRapidityMu", "P_{T} distribution vs Rapidity for #mu", nbinsRapidity,rapidityMin,rapidityMax, nbinsPt,ptMin,ptMax,fShouldSeparatePlusAndMinus);
206   
207   CreateTrackHisto(eventSelection,triggerClassName,centrality,"PtPhiMu", "P_{T} distribution vs phi for #mu", nbinsPhi,phiMin,phiMax, nbinsPt,ptMin,ptMax,fShouldSeparatePlusAndMinus);
208   
209   
210   CreateTrackHisto(eventSelection,triggerClassName,centrality,"PEtaMu", "P distribution for #mu",nbinsEta,etaMin,etaMax,nbinsP,pMin,pMax,fShouldSeparatePlusAndMinus);
211   
212   Double_t chi2min = 0;
213   Double_t chi2max = 20;
214   Int_t nbinchi2 = GetNbins(chi2min,chi2max,0.05);
215   
216   CreateTrackHisto(eventSelection, triggerClassName, centrality, "Chi2Mu", "chisquare per NDF #mu", nbinchi2, chi2min, chi2max,-1, 0.0, 0.0, fShouldSeparatePlusAndMinus);
217   
218   Double_t xmin = 0;
219   Double_t xmax = 150;
220   Int_t nbins = GetNbins(xmin,xmax,2.0);
221   
222   CreateTrackHisto(eventSelection,triggerClassName,centrality,"dcaP23Mu","#mu DCA vs P for 2-3 degrees;P (GeV);DCA (cm)",nbinsP,pMin,pMax,nbins,xmin,xmax,fShouldSeparatePlusAndMinus);
223   
224   CreateTrackHisto(eventSelection,triggerClassName,centrality,"dcaP310Mu","#mu DCA vs P for 3-10 degrees;P (GeV);DCA (cm)",nbinsP,pMin,pMax,nbins,xmin,xmax,fShouldSeparatePlusAndMinus);
225   
226   CreateTrackHisto(eventSelection,triggerClassName,centrality,"dcaPwPtCut23Mu","#mu DCA vs P for 2-3 degrees with Pt Cut;P (GeV);DCA (cm)",nbinsP,pMin,pMax,nbins,xmin,xmax,fShouldSeparatePlusAndMinus);
227   
228   CreateTrackHisto(eventSelection,triggerClassName,centrality,"dcaPwPtCut310Mu","#mu DCA vs P for 3-10 degrees with Pt Cut;P (GeV);DCA (cm)",nbinsP,pMin,pMax,nbins,xmin,xmax,fShouldSeparatePlusAndMinus);
229   
230 }
231
232 //_____________________________________________________________________________
233 void AliAnalysisMuMuSingle::FillHistosForTrack(const char* eventSelection,
234                                                const char* triggerClassName,
235                                                const char* centrality,
236                                                const char* trackCutName,
237                                                const AliVParticle& track)
238 {
239   /// Fill histograms for one track
240   
241   if (!AliAnalysisMuonUtility::IsMuonTrack(&track) ) return;
242   
243   if ( HasMC() )
244   {
245     MuonTrackCuts()->SetIsMC();
246   }
247   
248   TLorentzVector p(track.Px(),track.Py(),track.Pz(),
249                    TMath::Sqrt(AliAnalysisMuonUtility::MuonMass2()+track.P()*track.P()));
250   
251   
252   TString charge("");
253   
254   if ( ShouldSeparatePlusAndMinus() )
255   {
256     if ( track.Charge() < 0 )
257     {
258       charge = "Minus";
259     }
260     else
261     {
262       charge = "Plus";
263     }
264   }
265   
266   Double_t dca = EAGetTrackDCA(track);
267   
268   Double_t theta = AliAnalysisMuonUtility::GetThetaAbsDeg(&track);
269   
270   if (!IsHistogramDisabled("Chi2MatchTrigger"))
271   {
272     Histo(eventSelection,triggerClassName,centrality,trackCutName,"Chi2MatchTrigger")->Fill(AliAnalysisMuonUtility::GetChi2MatchTrigger(&track));
273   }
274   
275   if (!IsHistogramDisabled("EtaRapidityMu*"))
276   {
277     Histo(eventSelection,triggerClassName,centrality,trackCutName,Form("EtaRapidityMu%s",charge.Data()))->Fill(p.Rapidity(),p.Eta());
278   }
279   
280   if (!IsHistogramDisabled("PtEtaMu*"))
281   {
282     Histo(eventSelection,triggerClassName,centrality,trackCutName,Form("PtEtaMu%s",charge.Data()))->Fill(p.Eta(),p.Pt());
283   }
284   
285   if (!IsHistogramDisabled("PtRapidityMu*"))
286   {
287     Histo(eventSelection,triggerClassName,centrality,trackCutName,Form("PtRapidityMu%s",charge.Data()))->Fill(p.Rapidity(),p.Pt());
288   }
289   
290   if (!IsHistogramDisabled("PEtaMu*"))
291   {
292     Histo(eventSelection,triggerClassName,centrality,trackCutName,Form("PEtaMu%s",charge.Data()))->Fill(p.Eta(),p.P());
293   }
294   
295   if (!IsHistogramDisabled("PtPhiMu*"))
296   {
297     Histo(eventSelection,triggerClassName,centrality,trackCutName,Form("PtPhiMu%s",charge.Data()))->Fill(p.Phi(),p.Pt());
298   }
299   
300   if (!IsHistogramDisabled("Chi2Mu*"))
301   {
302     Histo(eventSelection,triggerClassName,centrality,trackCutName,Form("Chi2Mu%s",charge.Data()))->Fill(AliAnalysisMuonUtility::GetChi2perNDFtracker(&track));
303   }
304   
305   if ( theta >= 2.0 && theta < 3.0 )
306   {
307     
308     if (!IsHistogramDisabled("dcaP23Mu*"))
309     {
310       Histo(eventSelection,triggerClassName,centrality,trackCutName,Form("dcaP23Mu%s",charge.Data()))->Fill(p.P(),dca);
311     }
312     
313     if ( p.Pt() > 2 )
314     {
315       if (!IsHistogramDisabled("dcaPwPtCut23Mu*"))
316       {
317         Histo(eventSelection,triggerClassName,centrality,trackCutName,Form("dcaPwPtCut23Mu%s",charge.Data()))->Fill(p.P(),dca);
318       }
319     }
320   }
321   else if ( theta >= 3.0 && theta < 10.0 )
322   {
323     if (!IsHistogramDisabled("dcaP310Mu*"))
324     {
325       Histo(eventSelection,triggerClassName,centrality,trackCutName,Form("dcaP310Mu%s",charge.Data()))->Fill(p.P(),dca);
326     }
327     if ( p.Pt() > 2 )
328     {
329       if (!IsHistogramDisabled("dcaPwPtCut310Mu*"))
330       {
331         Histo(eventSelection,triggerClassName,centrality,trackCutName,Form("dcaPwPtCut310Mu%s",charge.Data()))->Fill(p.P(),dca);
332       }
333     }
334   }
335 }
336
337 //_____________________________________________________________________________
338 AliMuonTrackCuts* AliAnalysisMuMuSingle::MuonTrackCuts()
339 {
340   /// Get (and create the first time) our internal track cuts
341   if (!fMuonTrackCuts)
342   {
343     fMuonTrackCuts = new AliMuonTrackCuts;
344     
345     fMuonTrackCuts->SetAllowDefaultParams(kTRUE);
346     
347     fMuonTrackCuts->SetFilterMask(AliMuonTrackCuts::kMuEta |
348                                   AliMuonTrackCuts::kMuThetaAbs |
349                                   AliMuonTrackCuts::kMuPdca |
350                                   AliMuonTrackCuts::kMuMatchApt |
351                                   AliMuonTrackCuts::kMuMatchLpt |
352                                   AliMuonTrackCuts::kMuMatchHpt |
353                                   AliMuonTrackCuts::kMuTrackChiSquare);
354     
355   }
356   
357   return fMuonTrackCuts;
358 }
359
360 //_____________________________________________________________________________
361 void AliAnalysisMuMuSingle::SetMuonTrackCuts(const AliMuonTrackCuts& trackCuts)
362 {
363   /// Set our muontrackcuts from external source
364   delete fMuonTrackCuts;
365   fMuonTrackCuts = static_cast<AliMuonTrackCuts*>(trackCuts.Clone());
366 }