]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/TRD/AliTRDinfoGen.cxx
fix monitoring histo binning/filling
[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 <TH1S.h>
39 #include <TPad.h>
40 #include <TFile.h>
41 #include <TTree.h>
42 #include <TROOT.h>
43 #include <TChain.h>
44 #include <TParticle.h>
45
46 #include "AliLog.h"
47 #include "AliAnalysisManager.h"
48 #include "AliGeomManager.h"
49 #include "AliCDBManager.h"
50 #include "AliCDBEntry.h"
51 #include "AliCDBPath.h"
52 #include "AliESDEvent.h"
53 #include "AliMCEvent.h"
54 #include "AliESDInputHandler.h"
55 #include "AliMCEventHandler.h"
56
57 #include "AliESDfriend.h"
58 #include "AliESDfriendTrack.h"
59 #include "AliESDHeader.h"
60 #include "AliESDtrack.h"
61 #include "AliESDv0.h"
62 #include "AliESDtrackCuts.h"
63 #include "AliMCParticle.h"
64 #include "AliPID.h"
65 #include "AliStack.h"
66 #include "AliTrackReference.h"
67 #include "TTreeStream.h"
68
69 #include <cstdio>
70 #include <climits>
71 #include <cstring>
72 #include <iostream>
73 #include <memory>
74
75 #include "AliTRDReconstructor.h"
76 #include "AliTRDrecoParam.h"
77 #include "AliTRDcalibDB.h"
78 #include "AliTRDtrackerV1.h"
79 #include "AliTRDgeometry.h"
80 #include "AliTRDtrackV1.h"
81 #include "AliTRDseedV1.h"
82 #include "AliTRDcluster.h"
83 #include "AliTRDinfoGen.h"
84 #include "AliTRDpwg1Helper.h"
85 #include "info/AliTRDtrackInfo.h"
86 #include "info/AliTRDeventInfo.h"
87 #include "info/AliTRDv0Info.h"
88 #include "info/AliTRDeventCuts.h"
89
90 ClassImp(AliTRDinfoGen)
91
92 const Float_t AliTRDinfoGen::fgkITS = 100.; // to be checked
93 const Float_t AliTRDinfoGen::fgkTPC = 290.;
94 const Float_t AliTRDinfoGen::fgkTRD = 365.;
95
96 const Float_t AliTRDinfoGen::fgkEvVertexZ = 15.;
97 const Int_t   AliTRDinfoGen::fgkEvVertexN = 1;
98
99 const Float_t AliTRDinfoGen::fgkTrkDCAxy  = 3.;
100 const Float_t AliTRDinfoGen::fgkTrkDCAz   = 10.;
101 const Int_t   AliTRDinfoGen::fgkNclTPC    = 70;
102 const Float_t AliTRDinfoGen::fgkPt        = 0.2;
103 const Float_t AliTRDinfoGen::fgkEta       = 0.9;
104 AliTRDReconstructor* AliTRDinfoGen::fgReconstructor(NULL);
105 AliTRDgeometry* AliTRDinfoGen::fgGeo(NULL);
106 //____________________________________________________________________
107 AliTRDinfoGen::AliTRDinfoGen()
108   :AliAnalysisTaskSE()
109   ,fEvTrigger(NULL)
110   ,fESDev(NULL)
111   ,fMCev(NULL)
112   ,fEventCut(NULL)
113   ,fTrackCut(NULL)
114   ,fV0Cut(NULL)
115   ,fOCDB("local://$ALICE_ROOT/OCDB")
116   ,fTrackInfo(NULL)
117   ,fEventInfo(NULL)
118   ,fV0Info(NULL)
119   ,fTracksBarrel(NULL)
120   ,fTracksSA(NULL)
121   ,fTracksKink(NULL)
122   ,fV0List(NULL)
123   ,fContainer(NULL)
124   ,fDebugStream(NULL)
125 {
126   //
127   // Default constructor
128   //
129   SetNameTitle("TRDinfoGen", "MC-REC TRD-track list generator");
130 }
131
132 //____________________________________________________________________
133 AliTRDinfoGen::AliTRDinfoGen(char* name)
134   :AliAnalysisTaskSE(name)
135   ,fEvTrigger(NULL)
136   ,fESDev(NULL)
137   ,fMCev(NULL)
138   ,fEventCut(NULL)
139   ,fTrackCut(NULL)
140   ,fV0Cut(NULL)
141   ,fOCDB("local://$ALICE_ROOT/OCDB")
142   ,fTrackInfo(NULL)
143   ,fEventInfo(NULL)
144   ,fV0Info(NULL)
145   ,fTracksBarrel(NULL)
146   ,fTracksSA(NULL)
147   ,fTracksKink(NULL)
148   ,fV0List(NULL)
149   ,fContainer(NULL)
150   ,fDebugStream(NULL)
151 {
152   //
153   // Default constructor
154   //
155   SetTitle("MC-REC TRD-track list generator");
156   DefineOutput(AliTRDpwg1Helper::kTracksBarrel, TObjArray::Class());
157   DefineOutput(AliTRDpwg1Helper::kTracksSA, TObjArray::Class());
158   DefineOutput(AliTRDpwg1Helper::kTracksKink, TObjArray::Class());
159   DefineOutput(AliTRDpwg1Helper::kEventInfo, AliTRDeventInfo::Class());
160   DefineOutput(AliTRDpwg1Helper::kV0List, TObjArray::Class());
161   DefineOutput(AliTRDpwg1Helper::kMonitor, TObjArray::Class()); // histogram list
162 }
163
164 //____________________________________________________________________
165 AliTRDinfoGen::~AliTRDinfoGen()
166 {
167 // Destructor
168   if(fgGeo) delete fgGeo;
169   if(fgReconstructor) delete fgReconstructor;
170   if(fDebugStream) delete fDebugStream;
171   if(fEvTrigger) delete fEvTrigger;
172   if(fV0Cut) delete fV0Cut;
173   if(fTrackCut) delete fTrackCut;
174   if(fEventCut) delete fEventCut;
175   if(fTrackInfo) delete fTrackInfo; fTrackInfo = NULL;
176   if(fEventInfo) delete fEventInfo; fEventInfo = NULL;
177   if(fV0Info) delete fV0Info; fV0Info = NULL;
178   if(fTracksBarrel){ 
179     fTracksBarrel->Delete(); delete fTracksBarrel;
180     fTracksBarrel = NULL;
181   }
182   if(fTracksSA){ 
183     fTracksSA->Delete(); delete fTracksSA;
184     fTracksSA = NULL;
185   }
186   if(fTracksKink){ 
187     fTracksKink->Delete(); delete fTracksKink;
188     fTracksKink = NULL;
189   }
190   if(fV0List){ 
191     fV0List->Delete(); 
192     delete fV0List;
193     fV0List = NULL;
194   }
195   if(fContainer){ 
196     fContainer->Delete(); 
197     delete fContainer;
198     fContainer = NULL;
199   }
200 }
201
202 //____________________________________________________________________
203 Bool_t AliTRDinfoGen::GetRefFigure(Int_t)
204 {
205 // General graphs for PWG1/TRD train
206   if(!gPad){
207     AliWarning("Please provide a canvas to draw results.");
208     return kFALSE;
209   }
210   fContainer->At(0)->Draw("bar");
211   return kTRUE;
212 }
213
214 //____________________________________________________________________
215 void AliTRDinfoGen::UserCreateOutputObjects()
216 {       
217   //
218   // Create Output Containers (TObjectArray containing 1D histograms)
219   //
220  
221   fTrackInfo = new AliTRDtrackInfo();
222   fEventInfo = new AliTRDeventInfo();
223   fV0Info    = new AliTRDv0Info();
224   fTracksBarrel = new TObjArray(200); fTracksBarrel->SetOwner(kTRUE);
225   fTracksSA = new TObjArray(20); fTracksSA->SetOwner(kTRUE);
226   fTracksKink = new TObjArray(20); fTracksKink->SetOwner(kTRUE);
227   fV0List = new TObjArray(10); fV0List->SetOwner(kTRUE);
228
229   // define general monitor
230   fContainer = new TObjArray(1); fContainer->SetOwner(kTRUE);
231   TH1 *h=new TH1I("hStat", "Run statistics;Observable;Entries", Int_t(kNObjects), -0.5, Float_t(kNObjects)-0.5);
232   TAxis *ax(h->GetXaxis());
233   ax->SetBinLabel(Int_t(kTracksESD) + 1, "ESD");
234   ax->SetBinLabel(Int_t(kTracksMC) + 1, "MC");
235   ax->SetBinLabel(Int_t(kV0) + 1, "V0");
236   ax->SetBinLabel(Int_t(kTPC) + 1, "TPC");
237   ax->SetBinLabel(Int_t(kTRDin) + 1, "TRDin");
238   ax->SetBinLabel(Int_t(kTRDout) + 1, "TRDout");
239   ax->SetBinLabel(Int_t(kBarrel) + 1, "Barrel");
240   ax->SetBinLabel(Int_t(kBarrelMC) + 1, "BarrelMC");
241   ax->SetBinLabel(Int_t(kSA) + 1, "SA");
242   ax->SetBinLabel(Int_t(kSAMC) + 1, "SAMC");
243   ax->SetBinLabel(Int_t(kKink) + 1, "Kink");
244   ax->SetBinLabel(Int_t(kKinkMC) + 1, "KinkMC");
245   ax->SetBinLabel(Int_t(kBarrelFriend) + 1, "BFriend");
246   ax->SetBinLabel(Int_t(kSAFriend) + 1, "SFriend");
247   fContainer->AddAt(h, 0);
248   PostData(AliTRDpwg1Helper::kMonitor, fContainer);
249 }
250
251 //____________________________________________________________________
252 Bool_t AliTRDinfoGen::Load(const Char_t *file, const Char_t *dir, const Char_t *name)
253 {
254 // Load data from performance file
255
256   if(!TFile::Open(file)){
257     AliWarning(Form("Couldn't open file %s.", file));
258     return kFALSE;
259   }
260   if(dir){
261     if(!gFile->cd(dir)){
262       AliWarning(Form("Couldn't cd to %s in %s.", dir, file));
263       return kFALSE;
264     }
265   }
266   TObjArray *o(NULL);
267   const Char_t *tn=(name ? name : GetName());
268   if(!(o = (TObjArray*)gDirectory->Get(tn))){
269     AliWarning(Form("Missing histogram container %s.", tn));
270     return kFALSE;
271   }
272   fContainer = (TObjArray*)o->Clone(GetName());
273   gFile->Close();
274   return kTRUE;
275 }
276
277 //____________________________________________________________________
278 void AliTRDinfoGen::UserExec(Option_t *){
279   //
280   // Run the Analysis
281   //
282
283   fTracksBarrel->Delete();
284   fTracksSA->Delete();
285   fTracksKink->Delete();
286   fV0List->Delete();
287   fEventInfo->Delete("");
288
289   fESDev = dynamic_cast<AliESDEvent*>(InputEvent());
290   if(!fESDev){
291     AliError("Failed retrieving ESD event");
292     return;
293   }
294   // WARNING
295   // This part may conflict with other detectors !!
296   if(!IsInitOCDB()){ 
297     AliInfo("Initializing OCDB ...");
298     // prepare OCDB access
299     AliCDBManager* ocdb = AliCDBManager::Instance();
300     ocdb->SetDefaultStorage(fOCDB.Data());
301     ocdb->SetRun(fESDev->GetRunNumber());
302     // create geo manager
303     AliCDBEntry* obj(NULL);
304     if(!(obj = ocdb->Get(AliCDBPath("GRP", "Geometry", "Data")))){
305       AliError("GEOMETRY failed initialization.");
306     } else {
307       AliGeomManager::SetGeometry((TGeoManager*)obj->GetObject());
308       AliGeomManager::GetNalignable("TRD");
309       AliGeomManager::ApplyAlignObjsFromCDB("TRD");
310     }
311     fgGeo = new AliTRDgeometry;
312     fgGeo->CreateClusterMatrixArray();
313     //init magnetic field
314     if(!TGeoGlobalMagField::Instance()->IsLocked() &&
315        !fESDev->InitMagneticField()){
316       AliError("MAGNETIC FIELD failed initialization.");
317     }
318
319     // set no of time bins
320     AliTRDtrackerV1::SetNTimeBins(AliTRDcalibDB::Instance()->GetNumberOfTimeBinsDCS());
321     AliInfo(Form("OCDB :  Loc[%s] Run[%d] TB[%d]", fOCDB.Data(), ocdb->GetRun(), AliTRDtrackerV1::GetNTimeBins()));
322
323     // set reco param valid for this run/event
324     AliInfo(Form("Initializing TRD reco params for EventSpecie[%d]...",
325       fESDev->GetEventSpecie()));
326     fgReconstructor = new AliTRDReconstructor();
327     if(!(obj = ocdb->Get(AliCDBPath("TRD", "Calib", "RecoParam")))){
328       AliError("RECO PARAM failed initialization.");
329       fgReconstructor->SetRecoParam(AliTRDrecoParam::GetLowFluxParam());
330     } else {
331       obj->PrintMetaData();
332       TObjArray *recos((TObjArray*)obj->GetObject());
333       for(Int_t ireco(0); ireco<recos->GetEntriesFast(); ireco++){
334         AliTRDrecoParam *reco((AliTRDrecoParam*)recos->At(ireco));
335         Int_t es(reco->GetEventSpecie());
336         if(!(es&fESDev->GetEventSpecie())) continue;
337         fgReconstructor->SetRecoParam(reco);
338         TString s;
339         if(es&AliRecoParam::kLowMult) s="LowMult";
340         else if(es&AliRecoParam::kHighMult) s="HighMult";
341         else if(es&AliRecoParam::kCosmic) s="Cosmic";
342         else if(es&AliRecoParam::kCalib) s="Calib";
343         else s="Unknown";
344         AliInfo(Form("Using reco params for %s", s.Data()));
345         break;
346       }
347     }
348     SetInitOCDB();
349   }
350
351   // link MC if available
352   fMCev = MCEvent();
353
354   // event selection : trigger cut
355   if(UseLocalEvSelection() && fEvTrigger){ 
356     Bool_t kTRIGGERED(kFALSE);
357     std::auto_ptr<TObjArray> trig(fEvTrigger->Tokenize(" "));
358     for(Int_t itrig=trig->GetEntriesFast(); itrig--;){
359       const Char_t *trigClass(((TObjString*)(*trig)[itrig])->GetName());
360       if(fESDev->IsTriggerClassFired(trigClass)) {
361         AliDebug(2, Form("Ev[%4d] Trigger[%s]", fESDev->GetEventNumberInFile(), trigClass));
362         kTRIGGERED = kTRUE;
363         break; 
364       }
365     }
366     if(!kTRIGGERED){ 
367       AliDebug(2, Form("Reject Ev[%4d] Trigger", fESDev->GetEventNumberInFile()));
368       return;
369     }
370     // select only physical events
371     if(fESDev->GetEventType() != 7){ 
372       AliDebug(2, Form("Reject Ev[%4d] EvType[%d]", fESDev->GetEventNumberInFile(), fESDev->GetEventType()));
373       return;
374     }
375   }
376
377   // if the required trigger is a collision trigger then apply event vertex cut
378   if(UseLocalEvSelection() && IsCollision()){
379     const AliESDVertex *vertex = fESDev->GetPrimaryVertex();
380     if(TMath::Abs(vertex->GetZv())<1.e-10 || 
381        TMath::Abs(vertex->GetZv())>fgkEvVertexZ || 
382        vertex->GetNContributors()<fgkEvVertexN) {
383       AliDebug(2, Form("Reject Ev[%4d] Vertex Zv[%f] Nv[%d]", fESDev->GetEventNumberInFile(), TMath::Abs(vertex->GetZv()), vertex->GetNContributors()));
384       return;
385     }
386   }
387
388   if(fEventCut && !fEventCut->IsSelected(fESDev, IsCollision())) return;
389
390   if(!fESDfriend){
391     AliError("Failed retrieving ESD friend event");
392     return;
393   }
394   if(HasMCdata() && !fMCev){
395     AliError("Failed retrieving MC event");
396     return;
397   }
398
399   new(fEventInfo)AliTRDeventInfo(fESDev->GetHeader(), const_cast<AliESDRun *>(fESDev->GetESDRun()));
400   
401   Bool_t *trackMap(NULL);
402   AliStack * mStack(NULL);
403   Int_t nTracksMC = HasMCdata() ? fMCev->GetNumberOfTracks() : 0, nTracksESD = fESDev->GetNumberOfTracks();
404   if(HasMCdata()){
405     mStack = fMCev->Stack();
406     if(!mStack){
407       AliError("Failed retrieving MC Stack");
408       return;
409     }
410     trackMap = new Bool_t[nTracksMC];
411     memset(trackMap, 0, sizeof(Bool_t) * nTracksMC);
412   }
413   
414   Double32_t dedx[100]; Int_t nSlices(0);
415   Int_t nTRDout(0), nTRDin(0), nTPC(0)
416        ,nclsTrklt
417        ,nBarrel(0), nSA(0), nKink(0)
418        ,nBarrelFriend(0), nSAFriend(0)
419        ,nBarrelMC(0), nSAMC(0), nKinkMC(0);
420   AliESDtrack *esdTrack = NULL;
421   AliESDfriendTrack *esdFriendTrack = NULL;
422   TObject *calObject = NULL;
423   AliTRDtrackV1 *track = NULL;
424   AliTRDseedV1 *tracklet = NULL;
425   AliTRDcluster *cl = NULL;
426
427
428   // LOOP 0 - over ESD v0s
429   Float_t bField(fESDev->GetMagneticField());
430   AliESDv0 *v0(NULL);
431   Int_t v0pid[AliPID::kSPECIES];
432   for(Int_t iv0(0); iv0<fESDev->GetNumberOfV0s(); iv0++){
433     if(!(v0 = fESDev->GetV0(iv0))) continue;
434     // register v0
435     if(fV0Cut) new(fV0Info) AliTRDv0Info(*fV0Cut);
436     else new(fV0Info) AliTRDv0Info();
437     fV0Info->SetMagField(bField);
438     fV0Info->SetV0tracks(fESDev->GetTrack(v0->GetPindex()), fESDev->GetTrack(v0->GetNindex()));
439     fV0Info->SetV0Info(v0);
440     fV0List->Add(new AliTRDv0Info(*fV0Info));//  kFOUND=kFALSE;
441   }
442
443
444   // LOOP 1 - over ESD tracks
445   AliTRDv0Info *v0info=NULL;
446   for(Int_t itrk = 0; itrk < nTracksESD; itrk++){
447     new(fTrackInfo) AliTRDtrackInfo();
448     esdTrack = fESDev->GetTrack(itrk);
449     AliDebug(3, Form("\n%3d ITS[%d] TPC[%d] TRD[%d]\n", itrk, esdTrack->GetNcls(0), esdTrack->GetNcls(1), esdTrack->GetNcls(2)));
450
451     if(esdTrack->GetStatus()&AliESDtrack::kTPCout) nTPC++;
452     if(esdTrack->GetStatus()&AliESDtrack::kTRDout) nTRDout++;
453     if(esdTrack->GetStatus()&AliESDtrack::kTRDin) nTRDin++;
454
455     // look at external track param
456     const AliExternalTrackParam *op = esdTrack->GetOuterParam();
457     Double_t xyz[3];
458     if(op){
459       op->GetXYZ(xyz);
460       op->Global2LocalPosition(xyz, op->GetAlpha());
461       AliDebug(3, Form("op @ X[%7.3f]\n", xyz[0]));
462     }
463
464     // read MC info
465     Int_t fPdg = -1;
466     Int_t label = -1; UInt_t alab=UINT_MAX;
467     if(HasMCdata()){
468       label = esdTrack->GetLabel(); 
469       alab = TMath::Abs(label);
470       // register the track
471       if(alab < UInt_t(nTracksMC)){ 
472         trackMap[alab] = kTRUE; 
473       } else { 
474         AliError(Form("MC label[%d] outside scope for Ev[%d] Trk[%d].", label, (Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), itrk));
475         continue; 
476       }
477       AliMCParticle *mcParticle = NULL; 
478       if(!(mcParticle = (AliMCParticle*) fMCev->GetTrack(alab))){
479         AliError(Form("MC particle label[%d] missing for Ev[%d] Trk[%d].", label, (Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), itrk));
480         continue;
481       }
482       fPdg = mcParticle->Particle()->GetPdgCode();
483       Int_t nRefs = mcParticle->GetNumberOfTrackReferences();
484       Int_t iref = 0; AliTrackReference *ref = NULL; 
485       while(iref<nRefs){
486         ref = mcParticle->GetTrackReference(iref);
487         if(ref->LocalX() > fgkTPC) break;
488         iref++;
489       }
490
491       fTrackInfo->SetMC();
492       fTrackInfo->SetPDG(fPdg);
493       fTrackInfo->SetPrimary(mcParticle->Particle()->IsPrimary());
494       fTrackInfo->SetLabel(label);
495       Int_t jref = iref;//, kref = 0;
496       while(jref<nRefs){
497         ref = mcParticle->GetTrackReference(jref);
498         if(ref->LocalX() > fgkTRD) break;
499         AliDebug(4, Form("  trackRef[%2d (%2d)] @ %7.3f OK", jref-iref, jref, ref->LocalX()));
500         fTrackInfo->AddTrackRef(ref);
501         jref++;
502       }
503       AliDebug(3, Form("NtrackRefs[%d(%d)]", fTrackInfo->GetNTrackRefs(), nRefs));
504     }
505
506     // copy some relevant info to TRD track info
507     fTrackInfo->SetStatus(esdTrack->GetStatus());
508     fTrackInfo->SetTrackId(esdTrack->GetID());
509     Double_t p[AliPID::kSPECIES]; esdTrack->GetTRDpid(p);
510     fTrackInfo->SetESDpid(p);
511     fTrackInfo->SetESDpidQuality(esdTrack->GetTRDntrackletsPID());
512     if(!nSlices) nSlices = esdTrack->GetNumberOfTRDslices();
513     memset(dedx, 0, 100*sizeof(Double32_t));
514     Int_t in(0);
515     for(Int_t il=0; il<AliTRDgeometry::kNlayer; il++)
516       for(Int_t is=0; is<nSlices; is++) 
517         dedx[in++]=esdTrack->GetTRDslice(il, is);
518     for(Int_t il=0; il<AliTRDgeometry::kNlayer; il++) dedx[in++]=esdTrack->GetTRDmomentum(il);
519     fTrackInfo->SetSlices(in, dedx);
520     fTrackInfo->SetNumberOfClustersRefit(esdTrack->GetNcls(2));
521     // some other Informations which we may wish to store in order to find problematic cases
522     fTrackInfo->SetKinkIndex(esdTrack->GetKinkIndex(0));
523     fTrackInfo->SetTPCncls(static_cast<UShort_t>(esdTrack->GetNcls(1)));
524     nclsTrklt = 0;
525   
526     // set V0pid info
527     for(Int_t iv(0); iv<fV0List->GetEntriesFast(); iv++){
528       if(!(v0info = (AliTRDv0Info*)fV0List->At(iv))) continue;
529       if(!v0info->GetV0Daughter(1) && !v0info->GetV0Daughter(-1)) continue;
530       if(!v0info->HasTrack(fTrackInfo)) continue;
531       memset(v0pid, 0, AliPID::kSPECIES*sizeof(Int_t));
532       fTrackInfo->SetV0();
533       for(Int_t is=AliPID::kSPECIES; is--;){v0pid[is] = v0info->GetPID(is, fTrackInfo);}
534       fTrackInfo->SetV0pid(v0pid);
535       fTrackInfo->SetV0();
536       //const AliTRDtrackInfo::AliESDinfo *ei = fTrackInfo->GetESDinfo();
537       break;
538     }
539
540     // read REC info
541     esdFriendTrack = fESDfriend->GetTrack(itrk);
542     if(esdFriendTrack){
543       Int_t icalib = 0;
544       while((calObject = esdFriendTrack->GetCalibObject(icalib++))){
545         if(strcmp(calObject->IsA()->GetName(),"AliTRDtrackV1") != 0) continue; // Look for the TRDtrack
546         if(!(track = dynamic_cast<AliTRDtrackV1*>(calObject))) break;
547         AliDebug(4, Form("TRD track OK"));
548         // Set the clusters to unused
549         for(Int_t ipl = 0; ipl < AliTRDgeometry::kNlayer; ipl++){
550           if(!(tracklet = track->GetTracklet(ipl))) continue;
551           tracklet->ResetClusterIter();
552           while((cl = tracklet->NextCluster())) cl->Use(0);
553         }
554         fTrackInfo->SetTrack(track);
555         break;
556       }
557       AliDebug(3, Form("Ntracklets[%d]\n", fTrackInfo->GetNTracklets()));
558     } else AliDebug(3, "No ESD friends");
559     if(op) fTrackInfo->SetOuterParam(op);
560
561     if(DebugLevel() >= 1){
562       AliTRDtrackInfo info(*fTrackInfo);
563       (*DebugStream()) << "trackInfo"
564       << "TrackInfo.=" << &info
565       << "\n";
566       info.Delete("");
567     }
568
569     ULong_t status(esdTrack->GetStatus());
570     if((status&AliESDtrack::kTPCout)){
571       if(!esdTrack->GetKinkIndex(0)){ // Barrel  Track Selection
572         Bool_t selected(kTRUE);
573         if(UseLocalTrkSelection()){
574           if(esdTrack->Pt() < fgkPt){ 
575             AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] Pt[%5.2f]", fESDev->GetEventNumberInFile(), itrk, esdTrack->Pt()));
576             selected = kFALSE;
577           }
578           if(selected && TMath::Abs(esdTrack->Eta()) > fgkEta){
579             AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] Eta[%5.2f]", fESDev->GetEventNumberInFile(), itrk, TMath::Abs(esdTrack->Eta())));
580             selected = kFALSE;
581           }
582           if(selected && esdTrack->GetTPCNcls() < fgkNclTPC){ 
583             AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] NclTPC[%d]", fESDev->GetEventNumberInFile(), itrk, esdTrack->GetTPCNcls()));
584             selected = kFALSE;
585           }
586           Float_t par[2], cov[3];
587           esdTrack->GetImpactParameters(par, cov);
588           if(IsCollision()){ // cuts on DCA
589             if(selected && TMath::Abs(par[0]) > fgkTrkDCAxy){ 
590               AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] DCAxy[%f]", fESDev->GetEventNumberInFile(), itrk, TMath::Abs(par[0])));
591               selected = kFALSE;
592             }
593             if(selected && TMath::Abs(par[1]) > fgkTrkDCAz){ 
594               AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] DCAz[%f]", fESDev->GetEventNumberInFile(), itrk, TMath::Abs(par[1])));
595               selected = kFALSE;
596             }
597           } else if(selected && fMCev && !fMCev->IsPhysicalPrimary(alab)){;
598             AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] Primary", fESDev->GetEventNumberInFile(), itrk));
599             selected = kFALSE;
600           }
601         }
602         if(fTrackCut && !fTrackCut->IsSelected(esdTrack)) selected = kFALSE;
603         if(selected){ 
604           fTracksBarrel->Add(new AliTRDtrackInfo(*fTrackInfo));
605           nBarrel++;
606           if(fTrackInfo->GetTrack()) 
607             nBarrelFriend++;
608         }
609       } else {
610         fTracksKink->Add(new AliTRDtrackInfo(*fTrackInfo));
611         nKink++;
612       }
613     } else if((status&AliESDtrack::kTRDout) && !(status&AliESDtrack::kTRDin)){ 
614       fTracksSA->Add(new AliTRDtrackInfo(*fTrackInfo));
615       nSA++;
616       if(fTrackInfo->GetTrack()) 
617         nSAFriend++;
618     }
619     fTrackInfo->Delete("");
620   }
621
622   // LOOP 2 - over MC tracks which are passing TRD where the track is not reconstructed
623   if(HasMCdata()){
624     AliDebug(10, "Output of the MC track map:");
625     for(Int_t itk = 0; itk < nTracksMC;  itk++) AliDebug(10, Form("trackMap[%d] = %s", itk, trackMap[itk] == kTRUE ? "TRUE" : "kFALSE"));
626   
627     for(Int_t itk = 0; itk < nTracksMC; itk++){
628       if(trackMap[itk]) continue;
629       AliMCParticle *mcParticle =  (AliMCParticle*) fMCev->GetTrack(TMath::Abs(itk));
630       Int_t fPdg = mcParticle->Particle()->GetPdgCode();
631       Int_t nRefs = mcParticle->GetNumberOfTrackReferences();
632       Int_t iref = 0; AliTrackReference *ref = NULL; 
633       Int_t nRefsTRD = 0;
634       new(fTrackInfo) AliTRDtrackInfo();
635       fTrackInfo->SetMC();
636       fTrackInfo->SetPDG(fPdg);
637       while(iref<nRefs){ // count TRD TR
638         Bool_t kIN(kFALSE);
639         ref = mcParticle->GetTrackReference(iref);
640         if(ref->LocalX() > fgkTPC && ref->LocalX() < fgkTRD){
641           fTrackInfo->AddTrackRef(ref);
642           nRefsTRD++;kIN=kTRUE;
643         }
644         AliDebug(4, Form("  trackRef[%2d] @ x[%7.3f] %s", iref, ref->LocalX(), kIN?"IN":"OUT"));
645         iref++;
646       }
647       if(!nRefsTRD){
648         // In this stage we at least require 1 hit inside TRD. What will be done with this tracks is a task for the 
649         // analysis job
650         fTrackInfo->Delete("");
651         continue;
652       }
653       fTrackInfo->SetPrimary(mcParticle->Particle()->IsPrimary());
654       fTrackInfo->SetLabel(itk);
655       if(DebugLevel() >= 1){
656         AliTRDtrackInfo info(*fTrackInfo);
657         (*DebugStream()) << "trackInfo"
658         << "TrackInfo.=" << &info
659         << "\n";
660         info.Delete("");
661       }
662       AliDebug(3, Form("Add MC track @ label[%d] nTRDrefs[%d].", itk, nRefsTRD));
663       // check where the track starts
664       ref = mcParticle->GetTrackReference(0);
665       if(ref->LocalX() < fgkITS){ 
666         fTracksBarrel->Add(new AliTRDtrackInfo(*fTrackInfo));
667         nBarrelMC++;
668       } else if(ref->LocalX() < fgkTPC) {
669         fTracksKink->Add(new AliTRDtrackInfo(*fTrackInfo));
670         nKinkMC++;
671       } else if(nRefsTRD>6){
672         fTracksSA->Add(new AliTRDtrackInfo(*fTrackInfo));
673         nSAMC++;
674       }
675       fTrackInfo->Delete("");
676     }
677     delete[] trackMap;
678   }
679   AliDebug(1, Form(
680     "\nEv[%3d] Tracks: ESD[%d] MC[%d] V0[%d]\n"
681     "        TPCout[%d] TRDin[%d] TRDout[%d]\n"
682     "        Barrel[%3d+%3d=%3d] SA[%2d+%2d=%2d] Kink[%2d+%2d=%2d]"
683     ,(Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), nTracksESD, nTracksMC, fV0List->GetEntries()
684     , nTPC, nTRDin, nTRDout
685     ,nBarrel, nBarrelMC, fTracksBarrel->GetEntries()
686     ,nSA, nSAMC, fTracksSA->GetEntries()
687     ,nKink, nKinkMC, fTracksKink->GetEntries()
688   ));
689   // save track statistics
690   TH1 *h((TH1S*)fContainer->At(0));
691   h->Fill(Float_t(kTracksESD), nTracksESD);
692   h->Fill(Float_t(kTracksMC), nTracksMC);
693   h->Fill(Float_t(kV0), fV0List->GetEntries());
694   h->Fill(Float_t(kTPC), nTPC);
695   h->Fill(Float_t(kTRDin), nTRDin);
696   h->Fill(Float_t(kTRDout), nTRDout);
697   h->Fill(Float_t(kBarrel), nBarrel);
698   h->Fill(Float_t(kBarrelMC), nBarrelMC);
699   h->Fill(Float_t(kSA), nSA);
700   h->Fill(Float_t(kSAMC), nSAMC);
701   h->Fill(Float_t(kKink), nKink);
702   h->Fill(Float_t(kKinkMC), nKinkMC);
703   h->Fill(Float_t(kBarrelFriend), nBarrelFriend);
704   h->Fill(Float_t(kSAFriend), nSAFriend);
705
706   PostData(AliTRDpwg1Helper::kTracksBarrel, fTracksBarrel);
707   PostData(AliTRDpwg1Helper::kTracksSA, fTracksSA);
708   PostData(AliTRDpwg1Helper::kTracksKink, fTracksKink);
709   PostData(AliTRDpwg1Helper::kEventInfo, fEventInfo);
710   PostData(AliTRDpwg1Helper::kV0List, fV0List);
711 }
712
713 //____________________________________________________________________
714 void AliTRDinfoGen::SetLocalV0Selection(const AliTRDv0Info *v0)
715 {
716 // Set V0 cuts from outside
717
718   if(!fV0Cut) fV0Cut = new AliTRDv0Info(*v0);
719   else new(fV0Cut) AliTRDv0Info(*v0);
720 }
721
722 //____________________________________________________________________
723 void AliTRDinfoGen::SetTrigger(const Char_t *trigger)
724 {
725   if(!fEvTrigger) fEvTrigger = new TString(trigger);
726   else (*fEvTrigger) = trigger;
727 }
728
729 //____________________________________________________________________
730 TTreeSRedirector* AliTRDinfoGen::DebugStream()
731 {
732 // Manage debug stream for task
733   if(!fDebugStream){
734     TDirectory *savedir = gDirectory;
735     fDebugStream = new TTreeSRedirector("TRD.DebugInfoGen.root");
736     savedir->cd();
737   }
738   return fDebugStream;
739 }
740
741
742