]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/TRD/AliTRDinfoGen.cxx
restructured the class AliTRDv0Info
[u/mrichter/AliRoot.git] / PWG1 / TRD / AliTRDinfoGen.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: AliTRDinfoGen.cxx 27496 2008-07-22 08:35:45Z cblume $ */
17
18 ////////////////////////////////////////////////////////////////////////////
19 //
20 //  Tender wagon for TRD performance/calibration train
21 //
22 // In this wagon the information from
23 //   - ESD
24 //   - Friends [if available]
25 //   - MC [if available]
26 // are grouped into AliTRDtrackInfo objects and fed to worker tasks
27 //
28 //  Authors:
29 //    Markus Fasel <M.Fasel@gsi.de>
30 //    Alexandru Bercuci <A.Bercuci@gsi.de>
31 //
32 ////////////////////////////////////////////////////////////////////////////
33
34 #include <TClonesArray.h>
35 #include <TObjArray.h>
36 #include <TObject.h>
37 #include <TString.h>
38 #include <TH1F.h>
39 #include <TFile.h>
40 #include <TTree.h>
41 #include <TROOT.h>
42 #include <TChain.h>
43 #include <TParticle.h>
44
45 #include "AliLog.h"
46 #include "AliAnalysisManager.h"
47 #include "AliESDEvent.h"
48 #include "AliMCEvent.h"
49 #include "AliESDInputHandler.h"
50 #include "AliMCEventHandler.h"
51
52 #include "AliESDfriend.h"
53 #include "AliESDfriendTrack.h"
54 #include "AliESDHeader.h"
55 #include "AliESDtrack.h"
56 #include "AliESDv0.h"
57 #include "AliESDtrackCuts.h"
58 #include "AliMCParticle.h"
59 #include "AliPID.h"
60 #include "AliStack.h"
61 #include "AliTRDtrackV1.h"
62 #include "AliTrackReference.h"
63 #include "AliTRDgeometry.h"
64 #include "AliTRDcluster.h"
65 #include "AliTRDseedV1.h"
66 #include "TTreeStream.h"
67
68 #include <cstdio>
69 #include <climits>
70 #include <cstring>
71 #include <iostream>
72
73 #include "AliTRDinfoGen.h"
74 #include "info/AliTRDtrackInfo.h"
75 #include "info/AliTRDeventInfo.h"
76 #include "info/AliTRDv0Info.h"
77 #include "info/AliTRDeventCuts.h"
78 #include "macros/AliTRDperformanceTrain.h"
79
80 ClassImp(AliTRDinfoGen)
81
82 const Float_t AliTRDinfoGen::fgkITS = 100.; // to be checked
83 const Float_t AliTRDinfoGen::fgkTPC = 290.;
84 const Float_t AliTRDinfoGen::fgkTRD = 365.;
85
86 const Float_t AliTRDinfoGen::fgkEvVertexZ = 15.;
87 const Int_t   AliTRDinfoGen::fgkEvVertexN = 1;
88
89 const Float_t AliTRDinfoGen::fgkTrkDCAxy  = 40.;
90 const Float_t AliTRDinfoGen::fgkTrkDCAz   = 15.;
91 const Int_t   AliTRDinfoGen::fgkNclTPC    = 100;
92 const Float_t AliTRDinfoGen::fgkPt        = 0.2;
93 const Float_t AliTRDinfoGen::fgkEta       = 0.9;
94
95 //____________________________________________________________________
96 AliTRDinfoGen::AliTRDinfoGen()
97   :AliAnalysisTaskSE()
98   ,fEvTrigger(NULL)
99   ,fESDev(NULL)
100   ,fMCev(NULL)
101   ,fEventCut(NULL)
102   ,fTrackCut(NULL)
103   ,fV0Cut(NULL)
104   ,fTrackInfo(NULL)
105   ,fEventInfo(NULL)
106   ,fV0Info(NULL)
107   ,fTracksBarrel(NULL)
108   ,fTracksSA(NULL)
109   ,fTracksKink(NULL)
110   ,fV0List(NULL)
111   ,fDebugStream(NULL)
112 {
113   //
114   // Default constructor
115   //
116   SetNameTitle("infoGen", "MC-REC TRD-track list generator");
117 }
118
119 //____________________________________________________________________
120 AliTRDinfoGen::AliTRDinfoGen(char* name)
121   :AliAnalysisTaskSE(name)
122   ,fEvTrigger(NULL)
123   ,fESDev(NULL)
124   ,fMCev(NULL)
125   ,fEventCut(NULL)
126   ,fTrackCut(NULL)
127   ,fV0Cut(NULL)
128   ,fTrackInfo(NULL)
129   ,fEventInfo(NULL)
130   ,fV0Info(NULL)
131   ,fTracksBarrel(NULL)
132   ,fTracksSA(NULL)
133   ,fTracksKink(NULL)
134   ,fV0List(NULL)
135   ,fDebugStream(NULL)
136 {
137   //
138   // Default constructor
139   //
140   SetTitle("MC-REC TRD-track list generator");
141   DefineOutput(kTracksBarrel, TObjArray::Class());
142   DefineOutput(kTracksSA, TObjArray::Class());
143   DefineOutput(kTracksKink, TObjArray::Class());
144   DefineOutput(kEventInfo, AliTRDeventInfo::Class());
145   DefineOutput(kV0List, TObjArray::Class());
146 }
147
148 //____________________________________________________________________
149 AliTRDinfoGen::~AliTRDinfoGen()
150 {
151 // Destructor
152   
153   if(fDebugStream) delete fDebugStream;
154   if(fEvTrigger) delete fEvTrigger;
155   if(fV0Cut) delete fV0Cut;
156   if(fTrackCut) delete fTrackCut;
157   if(fEventCut) delete fEventCut;
158   if(fTrackInfo) delete fTrackInfo; fTrackInfo = NULL;
159   if(fEventInfo) delete fEventInfo; fEventInfo = NULL;
160   if(fV0Info) delete fV0Info; fV0Info = NULL;
161   if(fTracksBarrel){ 
162     fTracksBarrel->Delete(); delete fTracksBarrel;
163     fTracksBarrel = NULL;
164   }
165   if(fTracksSA){ 
166     fTracksSA->Delete(); delete fTracksSA;
167     fTracksSA = NULL;
168   }
169   if(fTracksKink){ 
170     fTracksKink->Delete(); delete fTracksKink;
171     fTracksKink = NULL;
172   }
173   if(fV0List){ 
174     fV0List->Delete(); 
175     delete fV0List;
176     fV0List = NULL;
177   }
178 }
179
180 //____________________________________________________________________
181 void AliTRDinfoGen::UserCreateOutputObjects()
182 {       
183   //
184   // Create Output Containers (TObjectArray containing 1D histograms)
185   //
186  
187   fTrackInfo = new AliTRDtrackInfo();
188   fEventInfo = new AliTRDeventInfo();
189   fV0Info    = new AliTRDv0Info();
190   fTracksBarrel = new TObjArray(200); fTracksBarrel->SetOwner(kTRUE);
191   fTracksSA = new TObjArray(20); fTracksSA->SetOwner(kTRUE);
192   fTracksKink = new TObjArray(20); fTracksKink->SetOwner(kTRUE);
193   fV0List = new TObjArray(10); fV0List->SetOwner(kTRUE);
194 }
195
196 //____________________________________________________________________
197 void AliTRDinfoGen::UserExec(Option_t *){
198   //
199   // Run the Analysis
200   //
201
202   fESDev = dynamic_cast<AliESDEvent*>(InputEvent());
203   fMCev = MCEvent();
204
205   if(!fESDev){
206     AliError("Failed retrieving ESD event");
207     return;
208   }
209
210   // event selection : trigger cut
211   if(UseLocalEvSelection() && fEvTrigger){ 
212     Bool_t kTRIGGERED(kFALSE);
213     const TObjArray *trig = fEvTrigger->Tokenize(" ");
214     for(Int_t itrig=trig->GetEntriesFast(); itrig--;){
215       const Char_t *trigClass(((TObjString*)(*trig)[itrig])->GetName());
216       if(fESDev->IsTriggerClassFired(trigClass)) {
217         AliDebug(2, Form("Ev[%4d] Trigger[%s]", fESDev->GetEventNumberInFile(), trigClass));
218         kTRIGGERED = kTRUE;
219         break; 
220       }
221     }
222     if(!kTRIGGERED){ 
223       AliDebug(2, Form("Reject Ev[%4d] Trigger", fESDev->GetEventNumberInFile()));
224       return;
225     }
226     // select only physical events
227     if(fESDev->GetEventType() != 7){ 
228       AliDebug(2, Form("Reject Ev[%4d] EvType[%d]", fESDev->GetEventNumberInFile(), fESDev->GetEventType()));
229       return;
230     }
231   }
232
233   // if the required trigger is a collision trigger then apply event vertex cut
234   if(UseLocalEvSelection() && IsCollision()){
235     const AliESDVertex *vertex = fESDev->GetPrimaryVertex();
236     if(TMath::Abs(vertex->GetZv())<1.e-10 || 
237        TMath::Abs(vertex->GetZv())>fgkEvVertexZ || 
238        vertex->GetNContributors()<fgkEvVertexN) {
239       AliDebug(2, Form("Reject Ev[%4d] Vertex Zv[%f] Nv[%d]", fESDev->GetEventNumberInFile(), TMath::Abs(vertex->GetZv()), vertex->GetNContributors()));
240       return;
241     }
242   }
243
244   if(fEventCut && !fEventCut->IsSelected(fESDev, IsCollision())) return;
245
246   if(!fESDfriend){
247     AliError("Failed retrieving ESD friend event");
248     return;
249   }
250   if(HasMCdata() && !fMCev){
251     AliError("Failed retrieving MC event");
252     return;
253   }
254
255   fTracksBarrel->Delete();
256   fTracksSA->Delete();
257   fTracksKink->Delete();
258   fV0List->Delete();
259   fEventInfo->Delete("");
260   new(fEventInfo)AliTRDeventInfo(fESDev->GetHeader(), const_cast<AliESDRun *>(fESDev->GetESDRun()));
261   
262   Bool_t *trackMap(NULL);
263   AliStack * mStack(NULL);
264   Int_t nTracksMC = HasMCdata() ? fMCev->GetNumberOfTracks() : 0, nTracksESD = fESDev->GetNumberOfTracks();
265   if(HasMCdata()){
266     mStack = fMCev->Stack();
267     if(!mStack){
268       AliError("Failed retrieving MC Stack");
269       return;
270     }
271     trackMap = new Bool_t[nTracksMC];
272     memset(trackMap, 0, sizeof(Bool_t) * nTracksMC);
273   }
274   
275   Double32_t dedx[100]; Int_t nSlices(0);
276   Int_t nTRDout(0), nTRDin(0), nTPC(0)
277        ,nclsTrklt
278        ,nBarrel(0), nSA(0), nKink(0)
279        ,nBarrelMC(0), nSAMC(0), nKinkMC(0);
280   AliESDtrack *esdTrack = NULL;
281   AliESDfriendTrack *esdFriendTrack = NULL;
282   TObject *calObject = NULL;
283   AliTRDtrackV1 *track = NULL;
284   AliTRDseedV1 *tracklet = NULL;
285   AliTRDcluster *cl = NULL;
286   // LOOP 0 - over ESD tracks
287   for(Int_t itrk = 0; itrk < nTracksESD; itrk++){
288     new(fTrackInfo) AliTRDtrackInfo();
289     esdTrack = fESDev->GetTrack(itrk);
290     AliDebug(3, Form("\n%3d ITS[%d] TPC[%d] TRD[%d]\n", itrk, esdTrack->GetNcls(0), esdTrack->GetNcls(1), esdTrack->GetNcls(2)));
291
292     if(esdTrack->GetStatus()&AliESDtrack::kTPCout) nTPC++;
293     if(esdTrack->GetStatus()&AliESDtrack::kTRDout) nTRDout++;
294     if(esdTrack->GetStatus()&AliESDtrack::kTRDin) nTRDin++;
295
296     // look at external track param
297     const AliExternalTrackParam *op = esdTrack->GetOuterParam();
298     Double_t xyz[3];
299     if(op){
300       op->GetXYZ(xyz);
301       op->Global2LocalPosition(xyz, op->GetAlpha());
302       AliDebug(3, Form("op @ X[%7.3f]\n", xyz[0]));
303     }
304
305     // read MC info
306     Int_t fPdg = -1;
307     Int_t label = -1; UInt_t alab=UINT_MAX;
308     if(HasMCdata()){
309       label = esdTrack->GetLabel(); 
310       alab = TMath::Abs(label);
311       // register the track
312       if(alab < UInt_t(nTracksMC)){ 
313         trackMap[alab] = kTRUE; 
314       } else { 
315         AliError(Form("MC label[%d] outside scope for Ev[%d] Trk[%d].", label, (Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), itrk));
316         continue; 
317       }
318       AliMCParticle *mcParticle = NULL; 
319       if(!(mcParticle = (AliMCParticle*) fMCev->GetTrack(alab))){
320         AliError(Form("MC particle label[%d] missing for Ev[%d] Trk[%d].", label, (Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), itrk));
321         continue;
322       }
323       fPdg = mcParticle->Particle()->GetPdgCode();
324       Int_t nRefs = mcParticle->GetNumberOfTrackReferences();
325       Int_t iref = 0; AliTrackReference *ref = NULL; 
326       while(iref<nRefs){
327         ref = mcParticle->GetTrackReference(iref);
328         if(ref->LocalX() > fgkTPC) break;
329         iref++;
330       }
331
332       fTrackInfo->SetMC();
333       fTrackInfo->SetPDG(fPdg);
334       fTrackInfo->SetPrimary(mcParticle->Particle()->IsPrimary());
335       Int_t jref = iref;//, kref = 0;
336       while(jref<nRefs){
337         ref = mcParticle->GetTrackReference(jref);
338         if(ref->LocalX() > fgkTRD) break;
339         AliDebug(4, Form("  trackRef[%2d (%2d)] @ %7.3f OK", jref-iref, jref, ref->LocalX()));
340         fTrackInfo->AddTrackRef(ref);
341         jref++;
342       }
343       AliDebug(3, Form("NtrackRefs[%d(%d)]", fTrackInfo->GetNTrackRefs(), nRefs));
344     }
345
346     // copy some relevant info to TRD track info
347     fTrackInfo->SetStatus(esdTrack->GetStatus());
348     fTrackInfo->SetTrackId(esdTrack->GetID());
349     Double_t p[AliPID::kSPECIES]; esdTrack->GetTRDpid(p);
350     fTrackInfo->SetESDpid(p);
351     fTrackInfo->SetESDpidQuality(esdTrack->GetTRDntrackletsPID());
352     if(!nSlices) nSlices = esdTrack->GetNumberOfTRDslices();
353     memset(dedx, 0, 100*sizeof(Double32_t));
354     Int_t in(0);
355     for(Int_t il=0; il<AliTRDgeometry::kNlayer; il++)
356       for(Int_t is=0; is<nSlices; is++) 
357         dedx[in++]=esdTrack->GetTRDslice(il, is);
358     for(Int_t il=0; il<AliTRDgeometry::kNlayer; il++) dedx[in++]=esdTrack->GetTRDmomentum(il);
359     fTrackInfo->SetSlices(in, dedx);
360     fTrackInfo->SetLabel(label);
361     fTrackInfo->SetNumberOfClustersRefit(esdTrack->GetNcls(2));
362     // some other Informations which we may wish to store in order to find problematic cases
363     fTrackInfo->SetKinkIndex(esdTrack->GetKinkIndex(0));
364     fTrackInfo->SetTPCncls(static_cast<UShort_t>(esdTrack->GetNcls(1)));
365     nclsTrklt = 0;
366   
367
368     // read REC info
369     esdFriendTrack = fESDfriend->GetTrack(itrk);
370     if(esdFriendTrack){
371       Int_t icalib = 0;
372       while((calObject = esdFriendTrack->GetCalibObject(icalib++))){
373         if(strcmp(calObject->IsA()->GetName(),"AliTRDtrackV1") != 0) continue; // Look for the TRDtrack
374         if(!(track = dynamic_cast<AliTRDtrackV1*>(calObject))) break;
375         AliDebug(4, Form("TRD track OK"));
376         // Set the clusters to unused
377         for(Int_t ipl = 0; ipl < AliTRDgeometry::kNlayer; ipl++){
378           if(!(tracklet = track->GetTracklet(ipl))) continue;
379           tracklet->ResetClusterIter();
380           while((cl = tracklet->NextCluster())) cl->Use(0);
381         }
382         fTrackInfo->SetTrack(track);
383         break;
384       }
385       AliDebug(3, Form("Ntracklets[%d]\n", fTrackInfo->GetNTracklets()));
386     } else AliDebug(3, "No ESD friends");
387     if(op) fTrackInfo->SetOuterParam(op);
388
389     if(DebugLevel() >= 1){
390       AliTRDtrackInfo info(*fTrackInfo);
391       (*DebugStream()) << "trackInfo"
392       << "TrackInfo.=" << &info
393       << "\n";
394       info.Delete("");
395     }
396
397     ULong_t status(esdTrack->GetStatus());
398     if((status&AliESDtrack::kTPCout)){
399       if(!esdTrack->GetKinkIndex(0)){ // Barrel  Track Selection
400         Bool_t selected(kTRUE);
401         if(UseLocalTrkSelection()){
402           if(esdTrack->Pt() < fgkPt){ 
403             AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] Pt[%5.2f]", fESDev->GetEventNumberInFile(), itrk, esdTrack->Pt()));
404             selected = kFALSE;
405           }
406           if(selected && TMath::Abs(esdTrack->Eta()) > fgkEta){
407             AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] Eta[%5.2f]", fESDev->GetEventNumberInFile(), itrk, TMath::Abs(esdTrack->Eta())));
408             selected = kFALSE;
409           }
410           if(selected && esdTrack->GetTPCNcls() < fgkNclTPC){ 
411             AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] NclTPC[%d]", fESDev->GetEventNumberInFile(), itrk, esdTrack->GetTPCNcls()));
412             selected = kFALSE;
413           }
414           Float_t par[2], cov[3];
415           esdTrack->GetImpactParameters(par, cov);
416           if(IsCollision()){ // cuts on DCA
417             if(selected && TMath::Abs(par[0]) > fgkTrkDCAxy){ 
418               AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] DCAxy[%f]", fESDev->GetEventNumberInFile(), itrk, TMath::Abs(par[0])));
419               selected = kFALSE;
420             }
421             if(selected && TMath::Abs(par[1]) > fgkTrkDCAz){ 
422               AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] DCAz[%f]", fESDev->GetEventNumberInFile(), itrk, TMath::Abs(par[1])));
423               selected = kFALSE;
424             }
425           } else if(selected && fMCev && !fMCev->IsPhysicalPrimary(alab)){;
426             AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] Primary", fESDev->GetEventNumberInFile(), itrk));
427             selected = kFALSE;
428           }
429         }
430         if(fTrackCut && !fTrackCut->IsSelected(esdTrack)) selected = kFALSE;
431         if(selected){ 
432           fTracksBarrel->Add(new AliTRDtrackInfo(*fTrackInfo));
433           nBarrel++;
434         }
435       } else {
436         fTracksKink->Add(new AliTRDtrackInfo(*fTrackInfo));
437         nKink++;
438       }
439     } else if((status&AliESDtrack::kTRDout) && !(status&AliESDtrack::kTRDin)){ 
440       fTracksSA->Add(new AliTRDtrackInfo(*fTrackInfo));
441       nSA++;
442     }
443     fTrackInfo->Delete("");
444   }
445
446   // LOOP 1 - over ESD v0s
447   Bool_t kFOUND(kFALSE);
448   Float_t bField(fESDev->GetMagneticField());
449   AliESDv0 *v0(NULL);
450   for(Int_t iv0(0); iv0<fESDev->GetNumberOfV0s(); iv0++){
451     if(!(v0 = fESDev->GetV0(iv0))) continue;
452
453     if(fV0Cut) new(fV0Info) AliTRDv0Info(*fV0Cut);
454     else new(fV0Info) AliTRDv0Info();
455     fV0Info->SetMagField(bField);
456     fV0Info->SetV0tracks(fESDev->GetTrack(v0->GetPindex()), fESDev->GetTrack(v0->GetNindex()));
457     fV0Info->SetV0Info(v0);
458     // purge V0 irelevant for TRD
459     for(Int_t ib(0); ib<nBarrel; ib++){ 
460       if(fV0Info->HasTrack((AliTRDtrackInfo*)fTracksBarrel->At(ib))){
461         kFOUND=kTRUE;
462         break;
463       }
464     }
465     if(!kFOUND){
466       for(Int_t ik(0); ik<nKink; ik++){ 
467         if(fV0Info->HasTrack((AliTRDtrackInfo*)fTracksKink->At(ik))){
468           kFOUND=kTRUE;
469           break;
470         }
471       }
472     }
473     if(!kFOUND) continue;
474     fV0List->Add(new AliTRDv0Info(*fV0Info));
475   }
476
477   // LOOP 2 - over MC tracks which are passing TRD where the track is not reconstructed
478   if(HasMCdata()){
479     AliDebug(10, "Output of the MC track map:");
480     for(Int_t itk = 0; itk < nTracksMC;  itk++) AliDebug(10, Form("trackMap[%d] = %s", itk, trackMap[itk] == kTRUE ? "TRUE" : "kFALSE"));
481   
482     for(Int_t itk = 0; itk < nTracksMC; itk++){
483       if(trackMap[itk]) continue;
484       AliMCParticle *mcParticle =  (AliMCParticle*) fMCev->GetTrack(TMath::Abs(itk));
485       Int_t fPdg = mcParticle->Particle()->GetPdgCode();
486       Int_t nRefs = mcParticle->GetNumberOfTrackReferences();
487       Int_t iref = 0; AliTrackReference *ref = NULL; 
488       Int_t nRefsTRD = 0;
489       new(fTrackInfo) AliTRDtrackInfo();
490       fTrackInfo->SetMC();
491       fTrackInfo->SetPDG(fPdg);
492       while(iref<nRefs){ // count TRD TR
493         Bool_t kIN(kFALSE);
494         ref = mcParticle->GetTrackReference(iref);
495         if(ref->LocalX() > fgkTPC && ref->LocalX() < fgkTRD){
496           fTrackInfo->AddTrackRef(ref);
497           nRefsTRD++;kIN=kTRUE;
498         }
499         AliDebug(4, Form("  trackRef[%2d] @ x[%7.3f] %s", iref, ref->LocalX(), kIN?"IN":"OUT"));
500         iref++;
501       }
502       if(!nRefsTRD){
503         // In this stage we at least require 1 hit inside TRD. What will be done with this tracks is a task for the 
504         // analysis job
505         fTrackInfo->Delete("");
506         continue;
507       }
508       fTrackInfo->SetPrimary(mcParticle->Particle()->IsPrimary());
509       fTrackInfo->SetLabel(itk);
510       if(DebugLevel() >= 1){
511         AliTRDtrackInfo info(*fTrackInfo);
512         (*DebugStream()) << "trackInfo"
513         << "TrackInfo.=" << &info
514         << "\n";
515         info.Delete("");
516       }
517       AliDebug(3, Form("Add MC track @ label[%d] nTRDrefs[%d].", itk, nRefsTRD));
518       // check where the track starts
519       ref = mcParticle->GetTrackReference(0);
520       if(ref->LocalX() < fgkITS){ 
521         fTracksBarrel->Add(new AliTRDtrackInfo(*fTrackInfo));
522         nBarrelMC++;
523       } else if(ref->LocalX() < fgkTPC) {
524         fTracksKink->Add(new AliTRDtrackInfo(*fTrackInfo));
525         nKinkMC++;
526       } else if(nRefsTRD>6){
527         fTracksSA->Add(new AliTRDtrackInfo(*fTrackInfo));
528         nSAMC++;
529       }
530       fTrackInfo->Delete("");
531     }
532     delete[] trackMap;
533   }
534   AliDebug(2, Form(
535     "\nEv[%3d] Tracks: ESD[%d] MC[%d] V0[%d]\n"
536     "        TPCout[%d] TRDin[%d] TRDout[%d]\n"
537     "        Barrel[%3d+%3d=%3d] SA[%2d+%2d=%2d] Kink[%2d+%2d=%2d]"
538     ,(Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), nTracksESD, nTracksMC, fV0List->GetEntries()
539     , nTPC, nTRDin, nTRDout
540     ,nBarrel, nBarrelMC, fTracksBarrel->GetEntries()
541     ,nSA, nSAMC, fTracksSA->GetEntries()
542     ,nKink, nKinkMC, fTracksKink->GetEntries()
543   ));
544
545   PostData(kTracksBarrel, fTracksBarrel);
546   PostData(kTracksSA, fTracksSA);
547   PostData(kTracksKink, fTracksKink);
548   PostData(kEventInfo, fEventInfo);
549   PostData(kV0List, fV0List);
550 }
551
552 //____________________________________________________________________
553 void AliTRDinfoGen::SetLocalV0Selection(AliTRDv0Info *v0)
554 {
555 // Set V0 cuts from outside
556
557   if(!fV0Cut) fV0Cut = new AliTRDv0Info(*v0);
558   else new(fV0Cut) AliTRDv0Info(*v0);
559 }
560
561 //____________________________________________________________________
562 void AliTRDinfoGen::SetTrigger(const Char_t *trigger)
563 {
564   if(!fEvTrigger) fEvTrigger = new TString(trigger);
565   else (*fEvTrigger) = trigger;
566 }
567
568 //____________________________________________________________________
569 TTreeSRedirector* AliTRDinfoGen::DebugStream()
570 {
571   if(!fDebugStream){
572     TDirectory *savedir = gDirectory;
573     fDebugStream = new TTreeSRedirector("TRD.DebugInfoGen.root");
574     savedir->cd();
575   }
576   return fDebugStream;
577 }
578
579