]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/qaRec/AliTRDefficiencyMC.cxx
fix ordering of tasks
[u/mrichter/AliRoot.git] / TRD / qaRec / AliTRDefficiencyMC.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 /* $Id: AliTRDefficiencyMC.cxx 27496 2008-07-22 08:35:45Z cblume $ */
17
18 ////////////////////////////////////////////////////////////////////////////
19 //                                                                        //
20 //  Reconstruction QA                                                     //
21 //                                                                        //
22 //  Authors:                                                              //
23 //    Markus Fasel <M.Fasel@gsi.de>                                       //
24 //                                                                        //
25 ////////////////////////////////////////////////////////////////////////////
26
27 #include <TObjArray.h>
28 #include <TClonesArray.h>
29 #include <TProfile.h>
30 #include <TMath.h>
31 #include <TCanvas.h>
32 #include "TTreeStream.h"
33
34 #include "AliMagF.h"
35 #include "AliPID.h"
36 #include "AliTracker.h"
37 #include "AliMathBase.h"
38 #include "AliTrackReference.h"
39 #include "AliAnalysisManager.h"
40
41 #include "AliTRDcluster.h"
42 #include "AliTRDseedV1.h"
43 #include "AliTRDtrackV1.h"
44 #include "AliTRDtrackerV1.h"
45 #include "Cal/AliTRDCalPID.h"
46 #include "info/AliTRDtrackInfo.h"
47 #include "AliTRDinfoGen.h"
48 #include "AliTRDefficiencyMC.h"
49
50 ClassImp(AliTRDefficiencyMC)
51
52 //_____________________________________________________________________________
53 AliTRDefficiencyMC::AliTRDefficiencyMC()
54   :AliTRDrecoTask("EfficiencyMC", "Combined Tracking Efficiency")
55 {
56   //
57   // Default constructor
58   //
59 }
60
61
62 //_____________________________________________________________________________
63 void AliTRDefficiencyMC::CreateOutputObjects(){
64   //
65   // Create output objects
66   //
67
68   OpenFile(0, "RECREATE");
69   fContainer = Histos();
70
71 }
72
73 //_____________________________________________________________________________
74 void AliTRDefficiencyMC::Exec(Option_t *){
75   //
76   // Execute the task:
77   //
78   // Loop over TrackInfos
79   // 1st: check if there is a TRDtrack
80   // 2nd: put conditions on the track:
81   //      - check if we did not register it before
82   //      - check if we have Track References for the track 
83   // 3rd: Register track:
84   //      - accepted if both conditions are fulfilled
85   //      - contamination if at least one condition is not fulfilled
86   // 4th: check Monte-Carlo Track wheter findable or not if there is no TRD track in this track info
87   // 5th: register MC track as rejected if findable and not jet registered
88   // Match the registers accepted and rejected and clean register rejected
89   // Fill the histograms
90   //
91   const Int_t kArraySize = 10000;     // Store indices of track references in an array
92   Int_t index_accepted[kArraySize], index_rejected[kArraySize], index_contamination[kArraySize];
93   memset(index_accepted, 0, sizeof(Int_t) * kArraySize);
94   memset(index_rejected, 0, sizeof(Int_t) * kArraySize);
95   memset(index_contamination, 0, sizeof(Int_t) * kArraySize);
96   Int_t naccepted = 0, nrejected = 0, ncontamination = 0;
97   Bool_t isContamination = kFALSE;
98   
99   Int_t nTrackInfos = fTracks->GetEntriesFast();
100   AliTRDtrackV1 *TRDtrack = 0x0;
101   AliTRDtrackInfo *trkInf = 0x0;
102   for(Int_t itinf = 0; itinf < nTrackInfos; itinf++){
103     trkInf = dynamic_cast<AliTRDtrackInfo *>(fTracks->UncheckedAt(itinf));
104     if(!trkInf) continue;
105     if(trkInf->GetTrack() || trkInf->GetNumberOfClustersRefit()){
106       isContamination = IsRegistered(trkInf,index_accepted,naccepted);
107       if(!trkInf->GetNTrackRefs()){
108         // We reject the track since the Monte Carlo Information is missing
109         printf("Error: Track Reference missing for Track %d\n", trkInf->GetLabel());
110         isContamination = kTRUE;
111         // Debugging
112         if(TRDtrack && fDebugLevel > 5) FillStreamTrackWOMC(trkInf);
113       } 
114       if(isContamination){
115         // reject kink (we count these only once)
116         if(trkInf->GetKinkIndex()) continue;
117         // Register track as contamination
118         index_contamination[ncontamination++]=itinf;
119         continue;
120       }
121       // Accept track
122       if(fDebugLevel > 3)printf("Accept track\n");
123       // Register track as accepted
124       index_accepted[naccepted++] = itinf;
125     }else{
126       if(IsFindable(trkInf)){
127         // register track as rejected if not already registered there
128         // Attention:
129         // These track infos are not!!! registered as contamination
130         if(!IsRegistered(trkInf, index_rejected, nrejected)) index_rejected[nrejected++] = itinf;
131       }
132     }
133   }
134   // we have to check if the rejected tracks are registered as found
135   // a possible source for this:
136   // first the track is not found by the barrel tracking but it is later found
137   // by the stand alone tracking, then two track info objects with the same 
138   // label would be created
139   // Attention:
140   // these tracks are not! registered as contamination
141   Int_t tmprejected[kArraySize]; Int_t nrej = nrejected;
142   memcpy(tmprejected, index_rejected, sizeof(Int_t) * nrejected);
143   nrejected = 0;
144   for(Int_t irej = 0; irej < nrej; irej++){
145     trkInf = dynamic_cast<AliTRDtrackInfo *>(fTracks->At(tmprejected[irej]));
146     if(!IsRegistered(trkInf,index_accepted,naccepted)) index_rejected[nrejected++] = tmprejected[irej];
147   }
148   // Fill Histograms
149   FillHistograms(naccepted, &index_accepted[0], kAccepted);
150   FillHistograms(nrejected, &index_rejected[0], kRejected);
151   FillHistograms(ncontamination, &index_contamination[0], kContamination);
152   Int_t nall = naccepted + nrejected;
153   //if(fDebugLevel>=1)
154   printf("%3d Tracks: MC[%3d] TRD[%3d | %5.2f%%] \n", (Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), nall, naccepted, nall ? 1.E2*Float_t(naccepted)/Float_t(nall) : 0.);
155   printf("%3d Tracks: ALL[%3d] Contamination[%3d | %5.2f%%] \n", (Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), nall + ncontamination, ncontamination, nall ? 1.E2*Float_t(ncontamination)/Float_t(nall + ncontamination) : 0.);
156
157   PostData(0, fContainer);
158 }
159
160 //_____________________________________________________________________________
161 void AliTRDefficiencyMC::Terminate(Option_t *)
162 {
163   //
164   // Termination
165   //
166
167   if(fDebugStream){ 
168     delete fDebugStream;
169     fDebugStream = 0x0;
170     fDebugLevel = 0;
171   }
172
173   fContainer = dynamic_cast<TObjArray*>(GetOutputData(0));
174   if (!fContainer) {
175     Printf("ERROR: list not available");
176     return;
177   }
178 }
179
180 //_____________________________________________________________________________
181 Bool_t AliTRDefficiencyMC::PostProcess()
182 {
183   //
184   // Post Process 
185   //
186   // Change the histogram style
187   // For species histograms apply the colors connected with the given particle species
188   //
189   TH1 *histo = dynamic_cast<TH1 *>(fContainer->At(kEfficiencyHistogram));
190   histo->SetMarkerStyle(22);
191   histo->SetMarkerColor(kBlue);
192   histo->GetXaxis()->SetTitle("p [GeV/c]");
193   histo->GetXaxis()->SetMoreLogLabels();
194   histo->GetYaxis()->SetTitle("Efficiency [%]");
195   histo->GetYaxis()->SetRangeUser(0.99, 1.005);
196
197   histo = dynamic_cast<TH1 *>(fContainer->At(kContaminationHistogram));
198   histo->SetMarkerStyle(22);
199   histo->SetMarkerColor(kBlue);
200   histo->GetXaxis()->SetTitle("p [GeV/c]");
201   histo->GetXaxis()->SetMoreLogLabels();
202   histo->GetYaxis()->SetTitle("Contamination [%]");
203   
204   // Species Efficiency Histograms
205   for(Int_t ispec = 0; ispec < AliPID::kSPECIES; ispec++){
206     histo = dynamic_cast<TH1 *>(fContainer->At(kEfficiencySpeciesHistogram  + ispec));
207     histo->SetMarkerStyle(22);
208     histo->SetLineColor(AliTRDCalPID::GetPartColor(ispec));
209     histo->SetMarkerColor(AliTRDCalPID::GetPartColor(ispec));
210     histo->GetXaxis()->SetTitle("p [GeV/c]");
211     histo->GetXaxis()->SetMoreLogLabels();
212     histo->GetYaxis()->SetTitle("Efficiency [%]");
213     histo->GetYaxis()->SetRangeUser(0.99, 1.005);
214   }
215   
216   // Species Contamination Histograms
217   for(Int_t ispec = 0; ispec < AliPID::kSPECIES; ispec++){
218     histo = dynamic_cast<TH1 *>(fContainer->At(kContaminationSpeciesHistogram  + ispec));
219     histo->SetMarkerStyle(22);
220     histo->SetLineColor(AliTRDCalPID::GetPartColor(ispec));
221     histo->SetMarkerColor(AliTRDCalPID::GetPartColor(ispec));
222     histo->GetXaxis()->SetTitle("p [GeV/c]");
223     histo->GetXaxis()->SetMoreLogLabels();
224     histo->GetYaxis()->SetTitle("Contamination [%]");
225   }
226   
227   fNRefFigures = 6;
228   return kTRUE;
229 }
230
231 //_____________________________________________________________________________
232 Bool_t AliTRDefficiencyMC::GetRefFigure(Int_t ifig){
233   //
234   // Plot the histograms
235   //
236   if(ifig >= fNRefFigures) return kFALSE;
237   if(ifig < 2){
238     (dynamic_cast<TH1 *>(fContainer->At(ifig)))->Draw("e1");
239     return kTRUE;
240   }
241   switch(ifig){
242   case 2:
243     (dynamic_cast<TH1 *>(fContainer->At(kEfficiencySpeciesHistogram)))->Draw("e1");
244     for(Int_t ispec = 1; ispec < AliPID::kSPECIES; ispec++)
245       (dynamic_cast<TH1 *>(fContainer->At(kEfficiencySpeciesHistogram + ispec)))->Draw("e1same");
246     break;
247   case 3:
248     (dynamic_cast<TH1 *>(fContainer->At(kContaminationSpeciesHistogram)))->Draw("e1");
249     for(Int_t ispec = 1; ispec < AliPID::kSPECIES; ispec++)
250       (dynamic_cast<TH1 *>(fContainer->At(kContaminationSpeciesHistogram + ispec)))->Draw("e1same");
251     break;
252   case 4:
253     (dynamic_cast<TH1 *>(fContainer->At(kEfficiencyNoPID)))->Draw("e1");
254     break;
255   case 5:
256     (dynamic_cast<TH1 *>(fContainer->At(kContaminationNoPID)))->Draw("e1");
257     break;
258   }
259   return kTRUE;
260 }
261
262 //_____________________________________________________________________________
263 TObjArray *AliTRDefficiencyMC::Histos(){
264   //
265   // Create the histograms
266   //
267   const Int_t nbins = 11;
268
269   if(fContainer) return fContainer;
270   Float_t xbins[nbins+1] = {.5, .7, .9, 1.3, 1.7, 2.4, 3.5, 4.5, 5.5, 7., 9., 11.};
271   
272   TString species[AliPID::kSPECIES] = {"Pions", "Muons", "Electrons", "Kaons", "Protons"};
273   TString species_short[AliPID::kSPECIES] = {"Pi", "Mu", "El", "Ka", "Pr"};
274   
275   fContainer = new TObjArray();
276   fContainer->AddAt(new TProfile("trEffComb", "Combined Tracking Efficiency", nbins, xbins), kEfficiencyHistogram);
277   fContainer->AddAt(new TProfile("trContComb", "Combined Tracking Contamination", nbins, xbins), kContaminationHistogram);
278   for(Int_t ispec = 0; ispec < AliPID::kSPECIES; ispec++)
279     fContainer->AddAt(new TProfile(Form("trEffComb%s", species_short[ispec].Data()), Form("Combined Tracking Efficiency %s", species[ispec].Data()), nbins, xbins), kEfficiencySpeciesHistogram + ispec);
280   for(Int_t ispec = 0; ispec < AliPID::kSPECIES; ispec++)
281     fContainer->AddAt(new TProfile(Form("trContComb%s", species_short[ispec].Data()), Form("Combined Tracking Contamination %s", species[ispec].Data()), nbins, xbins), kContaminationSpeciesHistogram + ispec);
282   fContainer->AddAt(new TProfile("trEffCombNoPID", "Combined Tracking Efficiency", nbins, xbins), kEfficiencyNoPID);
283   fContainer->AddAt(new TProfile("trContCombNoPID", "Combined Tracking Contamination", nbins, xbins), kContaminationNoPID);
284   return fContainer;
285 }
286
287 //_____________________________________________________________________________
288 Bool_t AliTRDefficiencyMC::IsFindable(AliTRDtrackInfo *trkInf){
289   //
290   // Apply Cuts on the Monte Carlo track references
291   // return whether track is findable or not
292   //
293   const Float_t kAlpha = 0.349065850;
294   
295   if(fDebugLevel>10) printf("Analysing Track References\n");
296   // Check if track is findable
297   Double_t mom = 0.;
298   Float_t xmin = 10000.0, xmax = 0.0; 
299   Float_t ymin = 0.0, ymax = 0.0;
300   Float_t zmin = 0.0, zmax = 0.0;
301   Float_t lastx = 0.0, x = 0.0;
302   Int_t nLayers = 0;
303   Int_t sector[20];
304   AliTrackReference *trackRef = 0x0;
305   for(Int_t itr = 0; itr < trkInf->GetNTrackRefs(); itr++){
306     trackRef = trkInf->GetTrackRef(itr);
307     if(fDebugLevel>10) printf("%d. x[%f], y[%f], z[%f]\n", itr, trackRef->LocalX(), trackRef->LocalY(), trackRef->Z());
308     x = trackRef->LocalX(); 
309         
310     // Be Sure that we are inside TRD
311     if(x < AliTRDinfoGen::xTPC || x > AliTRDinfoGen::xTOF) continue;    
312     sector[itr] = Int_t(trackRef->Alpha()/kAlpha);
313     if(x < xmin){
314       xmin = trackRef->LocalX();
315       ymin = trackRef->LocalY();
316       zmin = trackRef->Z();
317       mom = trackRef->P();
318     } else if(x > xmax){
319       xmax = trackRef->LocalX();
320       ymax = trackRef->LocalY();
321       zmax = trackRef->Z();
322     }
323     if(itr > 0){
324       Float_t dist = TMath::Abs(x - lastx);
325       if(fDebugLevel>10) printf("x = %f, lastx = %f, dist = %f\n", x, lastx, dist);
326       if(TMath::Abs(dist - 3.7) < 0.1) nLayers++;       // ref(i+1) has to be larger than ref(i)
327     }
328     lastx = x;
329   }
330
331   // Apply cuts
332   Bool_t findable = kTRUE;
333   if(trkInf->GetNTrackRefs() > 2 && xmax > xmin){
334     if(mom < 0.55) findable = kFALSE;                                                                   // momentum cut at 0.6
335     Double_t yangle = (ymax -ymin)/(xmax - xmin);
336     Double_t zangle = (zmax -zmin)/(xmax - xmin);
337     if(fDebugLevel>10){
338       printf("track: y-Angle = %f, z-Angle = %f\n", yangle, zangle);
339       printf("nLayers = %d\n", nLayers);
340     }
341     if(TMath::ATan(TMath::Abs(yangle)) > 45.) findable = kFALSE;
342     if(TMath::ATan(TMath::Abs(zangle)) > 45.) findable = kFALSE;
343     if(nLayers < 4) findable = kFALSE;
344     if(!trkInf->IsPrimary()) findable = kFALSE;
345     Bool_t samesec = kTRUE;
346     for(Int_t iref = 1; iref < trkInf->GetNTrackRefs(); iref++)
347       if(sector[iref] != sector[0]) samesec = kFALSE;
348     if(!samesec) findable = kFALSE;             // Discard all tracks which are passing more than one sector
349     if(fDebugLevel){
350       Double_t trackAngle = TMath::ATan(yangle);
351       Bool_t primary = trkInf->IsPrimary();
352       (*fDebugStream) << "NotFoundTrack"
353         << "Momentum="  << mom
354         << "trackAngle="<< trackAngle
355         << "NLayers="   << nLayers
356         << "Primary="   << primary
357         << "\n";
358     }
359   }
360   else
361     findable = kFALSE;
362   return findable;
363 }
364
365 //_____________________________________________________________________________
366 void AliTRDefficiencyMC::FillHistograms(Int_t nTracks, Int_t *indices, FillingMode_t mode){
367   //
368   // Fill Histograms in three different modes:
369   // 1st tracks which are found and accepted
370   // 2nd tracks which are not found and not already counted
371   // 3rd contaminating tracks: either double counts (kinks!) or tracks with no MC hit inside TRD
372   //
373   const Int_t pid[AliPID::kSPECIES] = {211,13,11,321,2212};
374   Double_t trkmom = 0.;   // the track momentum
375   Int_t trkpid = -1;      // particle species
376   AliTRDtrackInfo *trkInf = 0x0;
377   for(Int_t itk = 0; itk < nTracks; itk++){
378     trkInf = dynamic_cast<AliTRDtrackInfo *>(fTracks->At(indices[itk]));
379     if(fDebugLevel > 2)printf("Accepted MC track: %d\n", trkInf->GetLabel());
380     if(trkInf->GetNTrackRefs()){
381       // use Monte-Carlo Information for Momentum and PID
382       trkmom = trkInf->GetTrackRef(0)->P();
383       trkpid = trkInf->GetPDG();
384     }else{
385       // Use TPC Momentum
386       trkmom = trkInf->GetTrack()->P();
387     }
388     switch(mode){
389       case kAccepted:
390         (dynamic_cast<TProfile *>(fContainer->At(kEfficiencyHistogram)))->Fill(trkmom, 1);
391         (dynamic_cast<TProfile *>(fContainer->At(kContaminationHistogram)))->Fill(trkmom, 0);
392         break;
393       case kRejected:
394         (dynamic_cast<TProfile *>(fContainer->At(kEfficiencyHistogram)))->Fill(trkmom, 0);
395         (dynamic_cast<TProfile *>(fContainer->At(kContaminationHistogram)))->Fill(trkmom, 0);
396         break;
397       case kContamination:
398         (dynamic_cast<TProfile *>(fContainer->At(kContaminationHistogram)))->Fill(trkmom, 1);
399         break;
400     }
401     // Fill species histogram
402     Int_t part_spec = -1;
403     for(Int_t ispec = 0; ispec < AliPID::kSPECIES; ispec++){
404       if(trkpid == pid[ispec]) part_spec = ispec;
405     }
406     if(part_spec >= 0){
407       switch(mode){
408         case kAccepted:
409           (dynamic_cast<TProfile *>(fContainer->At(kEfficiencySpeciesHistogram + part_spec)))->Fill(trkmom, 1);
410           (dynamic_cast<TProfile *>(fContainer->At(kContaminationSpeciesHistogram + part_spec)))->Fill(trkmom, 0);
411           break;
412         case kRejected:
413           (dynamic_cast<TProfile *>(fContainer->At(kEfficiencySpeciesHistogram + part_spec)))->Fill(trkmom, 0);          (dynamic_cast<TProfile *>(fContainer->At(kContaminationSpeciesHistogram + part_spec)))->Fill(trkmom, 0);
414           break;
415         case kContamination:
416           (dynamic_cast<TProfile *>(fContainer->At(kContaminationSpeciesHistogram + part_spec)))->Fill(trkmom, 1);
417           break;
418       }
419     } else {
420       // The particle Type is not registered
421       (dynamic_cast<TProfile *>(fContainer->At(kEfficiencyNoPID)))->Fill(trkmom, 1);
422       (dynamic_cast<TProfile *>(fContainer->At(kContaminationNoPID)))->Fill(trkmom, 1);
423     }
424   }
425 }
426
427 //_____________________________________________________________________________
428 void AliTRDefficiencyMC::FillStreamTrackWOMC(AliTRDtrackInfo *trkInf){
429   // fill debug stream
430   // we want to know:
431   //  1. The event number
432   //  2. The track label
433   //  3. The TRD track label
434   //  4. The frequency of the TRD Label
435   //  5. Momentum from TPC (NO STAND ALONE TRACK)
436   //  6. TPC Phi angle
437   //  7. the TRD track
438   //  8. Monte Carlo PID
439   //  9. We check the Labels of the TRD track according to them we search the maching Monte-Carlo track.
440   //     From the matching Monte-Carlo track we store trackRefs, phi and momentum
441   // 10. We may also want to keep the kink index
442   Double_t mom = trkInf->GetESDinfo()->GetOuterParam()->P();
443   Int_t event = (Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry();
444   Int_t label = trkInf->GetLabel();
445   Int_t kinkIndex = trkInf->GetKinkIndex();
446   Int_t pdg = trkInf->GetPDG();
447   Double_t TPCphi = trkInf->GetESDinfo()->GetOuterParam()->Phi();
448   Int_t TRDlabels[180]; // Container for the cluster labels
449   Int_t sortlabels[360];        // Cluster Labels sorted according their occurancy
450   AliTRDseedV1 *tracklet = 0x0;
451   AliTRDcluster *c = 0x0;
452   Int_t nclusters = 0x0;
453   AliTRDtrackV1 *TRDtrack = trkInf->GetTrack();
454   for(Int_t il = 0; il < AliTRDgeometry::kNlayer; il++){
455     tracklet = TRDtrack->GetTracklet(il);
456     if(!tracklet) continue;
457     tracklet->ResetClusterIter();
458     c = 0x0;
459     while((c = tracklet->NextCluster())) TRDlabels[nclusters++] = c->GetLabel(0);
460   }
461   // Determine Label and Frequency
462   AliMathBase::Freq(nclusters, const_cast<const Int_t *>(&TRDlabels[0]), &sortlabels[0], kTRUE);
463   Int_t TRDLabel = sortlabels[0];
464   Int_t freqTRD = sortlabels[1];
465   // find the track info object matching to the TRD track
466   AliTRDtrackInfo *realtrack = 0;
467   TObjArrayIter rtiter(fTracks);
468   while((realtrack = (AliTRDtrackInfo *)rtiter())){
469     if(realtrack->GetLabel() != TRDLabel) continue;
470     break;
471   }
472   TClonesArray trackRefs("AliTrackReference");
473   Int_t realPdg = -1;
474   Double_t realP = 0.;
475   Double_t realPhi = 0.;
476   if(realtrack){
477     // pack the track references into the trackRefsContainer
478     for(Int_t iref = 0; iref < realtrack->GetNTrackRefs(); iref++){
479     new(trackRefs[iref])AliTrackReference(*(realtrack->GetTrackRef(iref)));
480     }
481     realPdg = realtrack->GetPDG();
482     if(realtrack->GetNTrackRefs()){
483       realP = realtrack->GetTrackRef(0)->P();
484       realPhi = realtrack->GetTrackRef(0)->Phi();
485     }
486   }
487   (*fDebugStream) << "TrackingEffMCfake"
488     << "Event=" << event
489     << "Label=" << label
490     << "TRDLabel=" << TRDLabel
491     << "FreqTRDlabel=" << freqTRD
492     << "TPCp="  << mom
493     << "TPCphi=" << TPCphi
494     << "TRDtrack=" << TRDtrack
495     << "PDG="   << pdg
496     << "TrackRefs=" << &trackRefs
497     << "RealPDG=" << realPdg
498     << "RealP=" << realP
499     << "RealPhi" << realPhi
500     << "KinkIndex=" << kinkIndex
501     << "\n";
502 }
503
504 //_____________________________________________________________________________
505 Bool_t AliTRDefficiencyMC::IsRegistered(AliTRDtrackInfo *trkInf, Int_t *indices, Int_t nTracks){
506   //
507   // Checks if track is registered in a given mode
508   //
509   Bool_t isRegistered = kFALSE;
510   for(Int_t il = 0; il < nTracks; il++){
511     if((dynamic_cast<AliTRDtrackInfo *>(fTracks->At(indices[il])))->GetLabel() == trkInf->GetLabel()){
512       isRegistered = kTRUE;        
513       break;
514     }
515   }
516   return isRegistered;
517 }
518