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