]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/Correlations/DPhi/FourierDecomposition/AliMuonEffMC.cxx
update DHC task (Constantin Loizides <cloizides@lbl.gov>)
[u/mrichter/AliRoot.git] / PWGCF / Correlations / DPhi / FourierDecomposition / AliMuonEffMC.cxx
1 // MUON track QA referring AliMuonEffMC.cxx
2 // Author : Saehanseul Oh
3
4 #include "AliMuonEffMC.h"
5
6 #include <TList.h>
7 #include <TH1D.h>
8 #include <TH2F.h>
9 #include <TH3F.h>
10 #include <THn.h>
11 #include <TChain.h>
12 #include <TFile.h>
13
14 #include "AliAnalysisManager.h"
15 #include "AliAnalysisTask.h"
16 #include "AliESDVertex.h"
17 #include "AliESDEvent.h"
18 #include "AliAODEvent.h"
19 #include "AliESDMuonTrack.h"
20 #include "AliAODTrack.h"
21 #include "AliESDVertex.h"
22 #include "AliAODVertex.h"
23 #include "AliCentrality.h"
24 #include "AliVParticle.h"
25 #include "AliMCParticle.h"
26 #include "AliAnalysisHelperJetTasks.h"
27 #include "AliMCEvent.h"
28
29 using std::cout;
30 using std::endl;
31
32 ClassImp(AliMuonEffMC)
33
34 //________________________________________________________________________
35 AliMuonEffMC::AliMuonEffMC() :
36   AliAnalysisTaskSE(), fESD(0), fAOD(0), fMC(0), fCentrality(99), fZVertex(99), fOutputList(0x0),      
37   fHEventStat(0), fHEvt(0x0), fIsMc(kTRUE), fMDProcess(kTRUE), fCentralityEstimator("V0M"),
38   fNEtaBins(15), fNpTBins(100), fNCentBins(5), fNZvtxBins(10), fNPhiBins(24),
39   fHMuonParGen(0x0), fHMuonDetGen(0x0), fHMuonDetRec(0x0), fHEtcDetRec(0x0)
40 {
41   // Constructor
42   DefineInput(0, TChain::Class());
43   DefineOutput(1, TList::Class());
44
45   for(Int_t i=0; i<4; i++)
46   {
47     fHMuMotherGenPt[i] = 0x0;
48     fHMuMotherRecPt[i] = 0x0;
49     fHMuMotherGenPhi[i] = 0x0;
50     fHMuMotherRecPhi[i] = 0x0;
51     fHMuMotherGenEta[i] = 0x0;
52     fHMuMotherRecEta[i] = 0x0;
53     fHMuDCA[i] = 0x0;
54   }
55 }
56
57 //________________________________________________________________________
58 AliMuonEffMC::AliMuonEffMC(const char *name) :
59   AliAnalysisTaskSE(name), fESD(0), fAOD(0), fMC(0), fCentrality(99), fZVertex(99), fOutputList(0x0),      
60   fHEventStat(0), fHEvt(0x0),  fIsMc(kTRUE), fMDProcess(kTRUE), fCentralityEstimator("V0M"),
61   fNEtaBins(15), fNpTBins(100), fNCentBins(5), fNZvtxBins(10), fNPhiBins(24),
62   fHMuonParGen(0x0), fHMuonDetGen(0x0), fHMuonDetRec(0x0), fHEtcDetRec(0x0)
63 {
64   // Constructor
65   DefineInput(0, TChain::Class());
66   DefineOutput(1, TList::Class());
67
68   for(Int_t i=0; i<4; i++)
69   {
70     fHMuMotherGenPt[i] = 0x0;
71     fHMuMotherRecPt[i] = 0x0;
72     fHMuMotherGenPhi[i] = 0x0;
73     fHMuMotherRecPhi[i] = 0x0;
74     fHMuMotherGenEta[i] = 0x0;
75     fHMuMotherRecEta[i] = 0x0;
76     fHMuDCA[i] = 0x0;
77   }
78 }
79
80 //________________________________________________________________________
81 AliMuonEffMC::~AliMuonEffMC()
82 {
83   //Destructor
84   if(fOutputList) delete fOutputList;
85 }
86
87 //________________________________________________________________________
88 void AliMuonEffMC::UserCreateOutputObjects()
89 {
90   // Create histograms
91   // Called once (per slave on PROOF!)
92
93   fOutputList = new TList();
94   fOutputList->SetOwner(1);
95  
96   fHEventStat = new TH1D("fHEventStat","Event statistics for analysis",18,0,18);
97   fHEventStat->GetXaxis()->SetBinLabel(1,"Event");
98   fHEventStat->GetXaxis()->SetBinLabel(2,"SelectedEvent");
99   fHEventStat->GetXaxis()->SetBinLabel(3,"File");
100   fHEventStat->GetXaxis()->SetBinLabel(4,"fSPHighpt");  //!Global Trigger Single plus High p_T
101   fHEventStat->GetXaxis()->SetBinLabel(5,"fSPAllpt");   //!Global Trigger Single plus All p_T
102   fHEventStat->GetXaxis()->SetBinLabel(6,"fSMLowpt");   //!Global Trigger Single minus Low p_T
103   fHEventStat->GetXaxis()->SetBinLabel(7,"fSMHighpt");  //!Global Trigger Single minus High p_T
104   fHEventStat->GetXaxis()->SetBinLabel(8,"fSMAllpt");   //!Global Trigger Single minus All p_T
105   fHEventStat->GetXaxis()->SetBinLabel(9,"fSULowpt");   //!Global Trigger Single undefined Low p_T
106   fHEventStat->GetXaxis()->SetBinLabel(10,"fSUHighpt"); //!Global Trigger Single undefined High p_T
107   fHEventStat->GetXaxis()->SetBinLabel(11,"fSUAllpt");  //!Global Trigger Single undefined All p_T
108   fHEventStat->GetXaxis()->SetBinLabel(12,"fUSLowpt");  //!Global Trigger UnlikeSign pair Low p_T
109   fHEventStat->GetXaxis()->SetBinLabel(13,"fUSHighpt"); //!Global Trigger UnlikeSign pair High p_T
110   fHEventStat->GetXaxis()->SetBinLabel(14,"fUSAllpt");  //!Global Trigger UnlikeSign pair All p_T
111   fHEventStat->GetXaxis()->SetBinLabel(15,"fLSLowpt");  //!Global Trigger LikeSign pair pair Low p_T
112   fHEventStat->GetXaxis()->SetBinLabel(16,"fLSHighpt"); //!Global Trigger LikeSign pair pair High p_T
113   fHEventStat->GetXaxis()->SetBinLabel(17,"fLSAllpt");  //!Global Trigger LikeSign pair pair All p_T
114   fHEventStat->GetXaxis()->SetBinLabel(18,"fSPLowpt");   //!Global Trigger Single plus Low p_T
115   fOutputList->Add(fHEventStat);
116
117   fHEvt = new TH2F("fHEvt", "Event-level variables; Zvtx; Cent", 30, -15, 15, 103, -2, 101);
118   fOutputList->Add(fHEvt);
119
120
121   Int_t iTrackBin[5];
122   Double_t* trackBins[5];
123   const char* trackAxisTitle[5];
124
125   // eta
126   Double_t etaBins[fNEtaBins+1];
127   for (Int_t i=0; i<=fNEtaBins; i++) { etaBins[i] = (Double_t)(-4.0 + 1.5/fNEtaBins*i); }
128   iTrackBin[0] = fNEtaBins;
129   trackBins[0] = etaBins;
130   trackAxisTitle[0] = "#eta";
131
132   // p_T
133   Double_t pTBins[fNpTBins+1];
134   for (Int_t i=0; i<=fNpTBins; i++) { pTBins[i] = (Double_t)(10.0/fNpTBins*i); }
135   iTrackBin[1] = fNpTBins;
136   trackBins[1] = pTBins;
137   trackAxisTitle[1] = "p_{T} (GeV/c)";
138
139   // centrality
140   Double_t CentBins[fNCentBins+1];
141   for (Int_t i=0; i<=fNCentBins; i++) { CentBins[i] = (Double_t)(100.0/fNCentBins*i); }
142   iTrackBin[2] = fNCentBins;
143   trackBins[2] = CentBins;
144   trackAxisTitle[2] = "Cent";
145
146   // Z-vertex
147   Double_t ZvtxBins[fNZvtxBins+1];
148   for (Int_t i=0; i<=fNZvtxBins; i++) { ZvtxBins[i] = (Double_t)(-10.0 + 20.0/fNZvtxBins*i); }
149   iTrackBin[3] = fNZvtxBins;
150   trackBins[3] = ZvtxBins;
151   trackAxisTitle[3] = "Zvtx";
152
153   // phi
154   Double_t phiBins[fNPhiBins+1];
155   for (Int_t i=0; i<=fNPhiBins; i++) { phiBins[i] = (Double_t)(TMath::TwoPi()/fNPhiBins * i); }
156   iTrackBin[4] = fNPhiBins;
157   trackBins[4] = phiBins;
158   trackAxisTitle[4] = "#phi";
159
160   // THn for tracking efficiency
161   fHMuonParGen = new THnF("fHMuonParGen", "", 5, iTrackBin, 0, 0);
162   for (Int_t i=0; i<5; i++)
163   {
164     fHMuonParGen->SetBinEdges(i, trackBins[i]);
165     fHMuonParGen->GetAxis(i)->SetTitle(trackAxisTitle[i]);
166   }
167   fHMuonParGen->Sumw2();
168   fOutputList->Add(fHMuonParGen);
169
170   fHMuonDetGen = (THnF*) fHMuonParGen->Clone("fHMuonDetGen");
171   fHMuonDetGen->Sumw2();
172   fOutputList->Add(fHMuonDetGen);
173
174   fHMuonDetRec = (THnF*) fHMuonParGen->Clone("fHMuonDetRec");
175   fHMuonDetRec->Sumw2();
176   fOutputList->Add(fHMuonDetRec);
177
178   fHEtcDetRec = (THnF*) fHMuonParGen->Clone("fHEtcDetRec");
179   fHEtcDetRec->Sumw2();
180   fOutputList->Add(fHEtcDetRec);
181
182   if(fMDProcess)
183   {
184     const char* MotherSpecies[4] = {"Pion","Kaon","D", "Etc"};
185     
186     for(Int_t i=0; i<4; i++)
187     {
188       fHMuMotherGenPt[i] = new TH2F(Form("fHMuMotherGenPt_%s",MotherSpecies[i]),";p_{T,muon}^{gen} (GeV/c);p_{T,mother}^{Truth} (GeV/c);",500, 0, 50, 500, 0, 50);
189       fOutputList->Add(fHMuMotherGenPt[i]);
190       
191       fHMuMotherRecPt[i] = new TH2F(Form("fHMuMotherRecPt_%s",MotherSpecies[i]),";p_{T,muon}^{rec} (GeV/c);p_{T,mother}^{Truth} (GeV/c);",500, 0, 50, 500, 0, 50);
192       fOutputList->Add(fHMuMotherRecPt[i]);
193       
194       fHMuMotherGenPhi[i] = new TH2F(Form("fHMuMotherGenPhi_%s",MotherSpecies[i]),";#phi_{gen};mother #phi;",100, 0, TMath::TwoPi(), 100, 0, TMath::TwoPi());
195       fOutputList->Add(fHMuMotherGenPhi[i]);
196       
197       fHMuMotherRecPhi[i] = new TH2F(Form("fHMuMotherRecPhi_%s",MotherSpecies[i]),";#phi_{rec};mother #phi;",100, 0, TMath::TwoPi(), 100, 0, TMath::TwoPi());
198       fOutputList->Add(fHMuMotherRecPhi[i]);
199       
200       fHMuMotherGenEta[i] = new TH2F(Form("fHMuMotherGenEta_%s",MotherSpecies[i]),";#eta_{gen};mother #eta;",100, -5., -1., 100, -5., -1.);
201       fOutputList->Add(fHMuMotherGenEta[i]);
202       
203       fHMuMotherRecEta[i] = new TH2F(Form("fHMuMotherRecEta_%s",MotherSpecies[i]),";#eta_{rec};mother #eta;",100, -5., -1., 100, -5., -1.);
204       fOutputList->Add(fHMuMotherRecEta[i]);
205       
206       fHMuDCA[i] =  new TH1F(Form("fHMuDCA_%s",MotherSpecies[i]), ";DCA", 100, 0, 50);
207       fOutputList->Add(fHMuDCA[i]);
208     }
209   }
210   PostData(1, fOutputList);
211 }
212
213 //________________________________________________________________________
214 void AliMuonEffMC::UserExec(Option_t *) 
215 {
216   // Main loop, Called for each event
217   Int_t ntrks = 0; // number of tracks in an event
218
219   if(((TString)InputEvent()->IsA()->GetName())=="AliAODEvent")
220   {
221     fAOD = dynamic_cast<AliAODEvent*>(InputEvent());
222     if (!fAOD) { AliError("AOD event not found. Nothing done!"); return; }
223     ntrks = fAOD->GetNTracks();
224   } 
225   else
226   {
227     fESD = dynamic_cast<AliESDEvent*>(InputEvent());
228     if (!fESD) { AliError("ESD event not found. Nothing done!"); return; }
229     ntrks = fESD->GetNumberOfMuonTracks();
230   }
231         
232   fHEventStat->Fill(0.5);
233   if(fIsMc) 
234   {
235     fMC = MCEvent();
236     if (!fMC) { AliError("MC event not avaliable."); return; }
237   }
238
239   // Centrality, vertex, other event variables...
240   if(fAOD)
241   { 
242     const AliAODVertex* vertex = fAOD->GetPrimaryVertex();
243     fZVertex = vertex->GetZ();
244     if(fAOD->GetCentrality())  fCentrality = fAOD->GetCentrality()->GetCentralityPercentile("fCentralityEstimator");
245   }
246   else if(fESD) 
247   {
248     const AliESDVertex* vertex = fESD->GetPrimaryVertex();
249     fZVertex = vertex->GetZ();
250     if(fESD->GetCentrality()) fCentrality = fESD->GetCentrality()->GetCentralityPercentile("fCentralityEstimator");
251   }   
252
253   if ((fESD && !VertexOk(fESD)) || (fAOD && !VertexOk(fAOD))) { AliInfo(Form("Event REJECTED. z = %.1f", fZVertex)); return; }
254   if (fCentrality > 100. || fCentrality < -1.5) { AliInfo(Form("Event REJECTED. fCentrality = %.1f", fCentrality)); return; }
255  
256   if(fCentrality < 0) fCentrality = 100.5;
257   // Fill Event histogram
258   fHEvt->Fill(fZVertex, fCentrality);
259   fHEventStat->Fill(1.5);
260  
261   ULong64_t trigword = 0;
262   if(fAOD) trigword=fAOD->GetTriggerMask();
263   else if(fESD) trigword=fESD->GetTriggerMask();
264  
265   if (trigword & 0x01) fHEventStat->Fill(17.5);
266   if (trigword & 0x02) fHEventStat->Fill(3.5);
267   if (trigword & 0x04) fHEventStat->Fill(4.5);
268   if (trigword & 0x08) fHEventStat->Fill(5.5);      
269   if (trigword & 0x010) fHEventStat->Fill(6.5);
270   if (trigword & 0x020) fHEventStat->Fill(7.5);
271   if (trigword & 0x040) fHEventStat->Fill(8.5);
272   if (trigword & 0x080) fHEventStat->Fill(9.5);
273   if (trigword & 0x100) fHEventStat->Fill(10.5);
274   if (trigword & 0x200) fHEventStat->Fill(11.5);
275   if (trigword & 0x400) fHEventStat->Fill(12.5);
276   if (trigword & 0x800) fHEventStat->Fill(13.5);
277   if (trigword & 0x1000) fHEventStat->Fill(14.5);
278   if (trigword & 0x2000) fHEventStat->Fill(15.5);
279   if (trigword & 0x4000) fHEventStat->Fill(16.5);
280
281   if(fIsMc)
282   {
283     for (Int_t iTrack = 0; iTrack<ntrks; iTrack++) 
284     { 
285       // loop on muon tracks
286       Int_t label = 0;
287       Double_t trackpt = 0;
288       Double_t tracketa = 0;
289       Double_t trackphi = 0;
290       Double_t dcavalue = 0;
291       if(fAOD)
292       {
293         AliAODTrack* muonTrack = (AliAODTrack*)fAOD->GetTrack(iTrack);
294         if(muonTrack)
295         {
296           if(!(IsGoodMUONtrack(*muonTrack)) || !(muonTrack->IsMuonTrack())) continue;
297           trackpt = muonTrack->Pt();
298           tracketa = muonTrack->Eta();
299           trackphi = muonTrack->Phi();
300           label =  TMath::Abs(muonTrack->GetLabel());
301         }
302       }
303       else if(fESD)
304       {
305         AliESDMuonTrack* muonTrack = fESD->GetMuonTrack(iTrack);
306         if(muonTrack)
307         {
308           if(!IsGoodMUONtrack(*muonTrack)) continue;
309           trackpt = muonTrack->Pt();
310           tracketa = muonTrack->Eta();
311           trackphi = muonTrack->Phi();
312           label =  TMath::Abs(muonTrack->GetLabel());
313           dcavalue = muonTrack->GetDCA();
314         }
315       }
316   
317
318       Double_t fillArrayDetRec[5] = { tracketa, trackpt, fCentrality, fZVertex, trackphi };
319       fHMuonDetRec->Fill(fillArrayDetRec); 
320         
321       AliMCParticle *McParticle  = (AliMCParticle*)fMC->GetTrack(label);
322       if(TMath::Abs(McParticle->PdgCode()) != 13) 
323       {
324         fHEtcDetRec->Fill(fillArrayDetRec); 
325         continue;
326       }
327       Double_t fillArrayDetGen[5] = { McParticle->Eta(), McParticle->Pt(), fCentrality, fZVertex, McParticle->Phi() };
328       fHMuonDetGen->Fill(fillArrayDetGen);
329         
330       if(fMDProcess)
331       {
332         if(McParticle->GetMother() > 0)
333         {
334           AliMCParticle *MotherParticle  = (AliMCParticle*)fMC->GetTrack(McParticle->GetMother());
335           Int_t motherlabel = TMath::Abs(MotherParticle->PdgCode());
336             
337           if(motherlabel==411 || motherlabel==413 || motherlabel==421 || motherlabel==423 || motherlabel==431 || motherlabel==433 || motherlabel==10413 || motherlabel==10411 || motherlabel==10423 || motherlabel==10421 || motherlabel==10433 || motherlabel==10431 || motherlabel==20413 || motherlabel==415 || motherlabel==20423 || motherlabel==425 || motherlabel==20433 || motherlabel==435)
338           {  
339             fHMuMotherGenPt[2]->Fill(McParticle->Pt(), MotherParticle->Pt());
340             fHMuMotherRecPt[2]->Fill(trackpt, MotherParticle->Pt());
341             fHMuMotherGenPhi[2]->Fill(McParticle->Phi(), MotherParticle->Phi());
342             fHMuMotherRecPhi[2]->Fill(trackphi, MotherParticle->Phi());
343             fHMuMotherGenEta[2]->Fill(McParticle->Eta(), MotherParticle->Eta());
344             fHMuMotherRecEta[2]->Fill(tracketa, MotherParticle->Eta());
345             if(fESD) fHMuDCA[2]->Fill(dcavalue);
346           }
347           else if(motherlabel==211) 
348             {
349               fHMuMotherGenPt[0]->Fill(McParticle->Pt(), MotherParticle->Pt());
350               fHMuMotherRecPt[0]->Fill(trackpt, MotherParticle->Pt());
351               fHMuMotherGenPhi[0]->Fill(McParticle->Phi(), MotherParticle->Phi());
352               fHMuMotherRecPhi[0]->Fill(trackphi, MotherParticle->Phi());
353               fHMuMotherGenEta[0]->Fill(McParticle->Eta(), MotherParticle->Eta());
354               fHMuMotherRecEta[0]->Fill(tracketa, MotherParticle->Eta());
355               if(fESD) fHMuDCA[0]->Fill(dcavalue);
356             }
357             else if(motherlabel==321)
358             {
359               fHMuMotherGenPt[1]->Fill(McParticle->Pt(), MotherParticle->Pt());
360               fHMuMotherRecPt[1]->Fill(trackpt, MotherParticle->Pt());
361               fHMuMotherGenPhi[1]->Fill(McParticle->Phi(), MotherParticle->Phi());
362               fHMuMotherRecPhi[1]->Fill(trackphi, MotherParticle->Phi());
363               fHMuMotherGenEta[1]->Fill(McParticle->Eta(), MotherParticle->Eta());
364               fHMuMotherRecEta[1]->Fill(tracketa, MotherParticle->Eta());
365               if(fESD) fHMuDCA[1]->Fill(dcavalue);
366             }   
367             else 
368             {
369               fHMuMotherGenPt[3]->Fill(McParticle->Pt(), MotherParticle->Pt());
370               fHMuMotherRecPt[3]->Fill(trackpt, MotherParticle->Pt());
371               fHMuMotherGenPhi[3]->Fill(McParticle->Phi(), MotherParticle->Phi());
372               fHMuMotherRecPhi[3]->Fill(trackphi, MotherParticle->Phi());
373               fHMuMotherGenEta[3]->Fill(McParticle->Eta(), MotherParticle->Eta());
374               fHMuMotherRecEta[3]->Fill(tracketa, MotherParticle->Eta());
375               if(fESD) fHMuDCA[3]->Fill(dcavalue);
376             }
377           }    
378         } // (mother hadron) : (daughter muon) QA 
379       } 
380     }
381     
382     for (Int_t ipart=0; ipart<fMC->GetNumberOfTracks(); ipart++)
383     {
384       AliMCParticle *McParticle  = (AliMCParticle*)fMC->GetTrack(ipart);
385       if(TMath::Abs(McParticle->PdgCode())==13) 
386       {
387         Double_t fillArrayParGen[5] = { McParticle->Eta(), McParticle->Pt(), fCentrality, fZVertex, McParticle->Phi() };
388         fHMuonParGen->Fill(fillArrayParGen);
389       }
390     }  
391
392   PostData(1, fOutputList);
393   return;
394 }
395
396 //________________________________________________________________________
397 void AliMuonEffMC::Terminate(Option_t *) 
398 {
399   // Draw result to the screen
400   // Called once at the end of the query
401
402 }
403
404 //________________________________________________________________________
405 Bool_t AliMuonEffMC::VertexOk(TObject* obj) const
406 {
407   // Modified from AliAnalyseLeadingTrackUE::VertexSelection()
408   
409   Int_t nContributors  = 0;
410   Double_t zVertex     = 999;
411   TString name("");
412   
413   if (obj->InheritsFrom("AliESDEvent")) {
414     AliESDEvent* esdevt = (AliESDEvent*) obj;
415     const AliESDVertex* vtx = esdevt->GetPrimaryVertex();
416     if (!vtx)
417       return 0;
418     nContributors = vtx->GetNContributors();
419     zVertex       = vtx->GetZ();
420     name          = vtx->GetName();
421   }
422   else if (obj->InheritsFrom("AliAODEvent")) { 
423     AliAODEvent* aodevt = (AliAODEvent*) obj;
424     if (aodevt->GetNumberOfVertices() < 1)
425       return 0;
426     const AliAODVertex* vtx = aodevt->GetPrimaryVertex();
427     nContributors = vtx->GetNContributors();
428     zVertex       = vtx->GetZ();
429     name          = vtx->GetName();
430   }
431   
432   // Reject if TPC-only vertex
433   if (name.CompareTo("TPCVertex")==0)
434     return kFALSE;
435   
436   // Check # contributors and range...
437   if( nContributors < 1 || TMath::Abs(zVertex) > 10 ) {
438     return kFALSE;
439   }
440   
441   return kTRUE;
442 }
443
444
445 //________________________________________________________________________
446 Bool_t AliMuonEffMC::IsGoodMUONtrack(AliESDMuonTrack &track)
447 {
448   // Applying track cuts for MUON tracks
449   if(!track.ContainTrackerData()) return kFALSE;
450
451   Double_t thetaTrackAbsEnd = TMath::ATan(track.GetRAtAbsorberEnd()/505.) * TMath::RadToDeg();
452   Double_t eta = track.Eta();
453
454   // Theta cut at absorber end
455   if(thetaTrackAbsEnd <= 2. || thetaTrackAbsEnd >= 10.) return kFALSE;
456   // Eta cut
457   if(eta <= -4. || eta >= -2.5) return kFALSE;
458   if(track.GetMatchTrigger() <= 0) return kFALSE;
459
460   return kTRUE;
461 }
462
463 //________________________________________________________________________
464 Bool_t AliMuonEffMC::IsGoodMUONtrack(AliAODTrack &track)
465 {
466   if (!track.IsMuonTrack()) return kFALSE;
467
468   Double_t dThetaAbs = TMath::ATan(track.GetRAtAbsorberEnd()/505.)
469                      * TMath::RadToDeg();
470   if ((dThetaAbs<2.) || (dThetaAbs>10.)) return kFALSE;
471
472   Double_t dEta = track.Eta();
473   if ((dEta<-4.) || (dEta>2.5)) return kFALSE;
474
475   if (track.GetMatchTrigger()<0.5) return kFALSE;
476
477   return kTRUE;
478 }