2dca0a9c3155ea1fe3a52b911ba48f3d972b1841
[u/mrichter/AliRoot.git] / PWG1 / TRD / AliTRDinfoGen.cxx
1 /**************************************************************************\r
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *\r
3 *                                                                        *\r
4 * Author: The ALICE Off-line Project.                                    *\r
5 * Contributors are mentioned in the code where appropriate.              *\r
6 *                                                                        *\r
7 * Permission to use, copy, modify and distribute this software and its   *\r
8 * documentation strictly for non-commercial purposes is hereby granted   *\r
9 * without fee, provided that the above copyright notice appears in all   *\r
10 * copies and that both the copyright notice and this permission notice   *\r
11 * appear in the supporting documentation. The authors make no claims     *\r
12 * about the suitability of this software for any purpose. It is          *\r
13 * provided "as is" without express or implied warranty.                  *\r
14 **************************************************************************/\r
15 \r
16 /* $Id: AliTRDinfoGen.cxx 27496 2008-07-22 08:35:45Z cblume $ */\r
17 \r
18 ////////////////////////////////////////////////////////////////////////////\r
19 //\r
20 //  Tender wagon for TRD performance/calibration train\r
21 //\r
22 // In this wagon the information from\r
23 //   - ESD\r
24 //   - Friends [if available]\r
25 //   - MC [if available]\r
26 // are grouped into AliTRDtrackInfo objects and fed to worker tasks\r
27 //\r
28 //  Authors:\r
29 //    Markus Fasel <M.Fasel@gsi.de>\r
30 //    Alexandru Bercuci <A.Bercuci@gsi.de>\r
31 //\r
32 ////////////////////////////////////////////////////////////////////////////\r
33 \r
34 #include <TClonesArray.h>\r
35 #include <TObjArray.h>\r
36 #include <TObject.h>\r
37 #include <TString.h>\r
38 #include <TVectorT.h>\r
39 #include <TH1.h>\r
40 #include <TH2.h>\r
41 #include <TCanvas.h>\r
42 #include <TPad.h>\r
43 #include <TFile.h>\r
44 #include <TTree.h>\r
45 #include <TROOT.h>\r
46 #include <TChain.h>\r
47 #include <TParticle.h>\r
48 \r
49 #include "AliLog.h"\r
50 #include "AliAnalysisManager.h"\r
51 #include "AliGeomManager.h"\r
52 #include "AliCDBManager.h"\r
53 #include "AliCDBStorage.h"\r
54 #include "AliCDBEntry.h"\r
55 #include "AliCDBPath.h"\r
56 #include "AliESDEvent.h"\r
57 #include "AliMCEvent.h"\r
58 #include "AliESDInputHandler.h"\r
59 #include "AliESDInputHandlerRP.h"\r
60 #include "AliMCEventHandler.h"\r
61 \r
62 #include "AliESDfriend.h"\r
63 #include "AliESDfriendTrack.h"\r
64 #include "AliESDHeader.h"\r
65 #include "AliESDRun.h"\r
66 #include "AliESDtrack.h"\r
67 #include "AliESDv0.h"\r
68 #include "AliESDtrackCuts.h"\r
69 #include "AliMCParticle.h"\r
70 #include "AliMultiplicity.h"\r
71 #include "AliCentrality.h"\r
72 #include "AliPID.h"\r
73 #include "AliStack.h"\r
74 #include "AliTrackReference.h"\r
75 #include "TTreeStream.h"\r
76 \r
77 #include <cstdio>\r
78 #include <climits>\r
79 #include <cstring>\r
80 #include <iostream>\r
81 #include <memory>\r
82 \r
83 #include "AliTRDReconstructor.h"\r
84 #include "AliTRDrecoParam.h"\r
85 #include "AliTRDcalibDB.h"\r
86 #include "AliTRDtrackerV1.h"\r
87 #include "AliTRDpadPlane.h"\r
88 #include "AliTRDgeometry.h"\r
89 #include "AliTRDtrackV1.h"\r
90 #include "AliTRDseedV1.h"\r
91 #include "AliTRDcluster.h"\r
92 #include "AliTRDinfoGen.h"\r
93 #include "AliTRDpwg1Helper.h"\r
94 #include "info/AliTRDtrackInfo.h"\r
95 #include "info/AliTRDeventInfo.h"\r
96 #include "info/AliTRDv0Info.h"\r
97 #include "info/AliTRDeventCuts.h"\r
98 \r
99 ClassImp(AliTRDinfoGen)\r
100 \r
101 const Float_t AliTRDinfoGen::fgkITS = 100.; // to be checked\r
102 const Float_t AliTRDinfoGen::fgkTPC = 290.;\r
103 const Float_t AliTRDinfoGen::fgkTRD = 365.;\r
104 \r
105 const Float_t AliTRDinfoGen::fgkTrkDCAxy  = 3.;\r
106 const Float_t AliTRDinfoGen::fgkTrkDCAz   = 10.;\r
107 const Int_t   AliTRDinfoGen::fgkNclTPC    = 70;\r
108 const Float_t AliTRDinfoGen::fgkPt        = 0.2;\r
109 const Float_t AliTRDinfoGen::fgkEta       = 0.9;\r
110 AliTRDReconstructor* AliTRDinfoGen::fgReconstructor(NULL);\r
111 AliTRDgeometry* AliTRDinfoGen::fgGeo(NULL);\r
112 //____________________________________________________________________\r
113 AliTRDinfoGen::AliTRDinfoGen()\r
114   :AliAnalysisTaskSE()\r
115   ,fESDev(NULL)\r
116   ,fMCev(NULL)\r
117   ,fEventCut(NULL)\r
118   ,fTrackCut(NULL)\r
119   ,fV0Cut(NULL)\r
120   ,fOCDB("local://$ALICE_ROOT/OCDB")\r
121   ,fTrackInfo(NULL)\r
122   ,fEventInfo(NULL)\r
123   ,fV0Info(NULL)\r
124   ,fTracksBarrel(NULL)\r
125   ,fTracksSA(NULL)\r
126   ,fTracksKink(NULL)\r
127   ,fV0List(NULL)\r
128   ,fClusters(NULL)\r
129   ,fContainer(NULL)\r
130   ,fRecos(NULL)\r
131   ,fDebugStream(NULL)\r
132 {\r
133   //\r
134   // Default constructor\r
135   //\r
136   SetNameTitle("TRDinfoGen", "MC-REC TRD-track list generator");\r
137 }\r
138 \r
139 //____________________________________________________________________\r
140 AliTRDinfoGen::AliTRDinfoGen(char* name)\r
141   :AliAnalysisTaskSE(name)\r
142   ,fESDev(NULL)\r
143   ,fMCev(NULL)\r
144   ,fEventCut(NULL)\r
145   ,fTrackCut(NULL)\r
146   ,fV0Cut(NULL)\r
147   ,fOCDB("local://$ALICE_ROOT/OCDB")\r
148   ,fTrackInfo(NULL)\r
149   ,fEventInfo(NULL)\r
150   ,fV0Info(NULL)\r
151   ,fTracksBarrel(NULL)\r
152   ,fTracksSA(NULL)\r
153   ,fTracksKink(NULL)\r
154   ,fV0List(NULL)\r
155   ,fClusters(NULL)\r
156   ,fContainer(NULL)\r
157   ,fRecos(NULL)\r
158   ,fDebugStream(NULL)\r
159 {\r
160   //\r
161   // Default constructor\r
162   //\r
163   SetTitle("MC-REC TRD-track list generator");\r
164   DefineOutput(AliTRDpwg1Helper::kTracksBarrel, TObjArray::Class());\r
165   DefineOutput(AliTRDpwg1Helper::kTracksSA,     TObjArray::Class());\r
166   DefineOutput(AliTRDpwg1Helper::kTracksKink,   TObjArray::Class());\r
167   DefineOutput(AliTRDpwg1Helper::kEventInfo,    AliTRDeventInfo::Class());\r
168   DefineOutput(AliTRDpwg1Helper::kV0List,       TObjArray::Class());\r
169   DefineOutput(AliTRDpwg1Helper::kClusters,     TObjArray::Class());\r
170   DefineOutput(AliTRDpwg1Helper::kMonitor,      TObjArray::Class()); // histogram list\r
171 }\r
172 \r
173 //____________________________________________________________________\r
174 AliTRDinfoGen::~AliTRDinfoGen()\r
175 {\r
176 // Destructor\r
177   if(fgGeo) delete fgGeo;\r
178   if(fgReconstructor) delete fgReconstructor;\r
179   if(fDebugStream) delete fDebugStream;\r
180   if(fV0Cut) delete fV0Cut;\r
181   if(fTrackCut) delete fTrackCut;\r
182   if(fEventCut) delete fEventCut;\r
183   if(fTrackInfo) delete fTrackInfo; fTrackInfo = NULL;\r
184   if(fEventInfo) delete fEventInfo; fEventInfo = NULL;\r
185   if(fV0Info) delete fV0Info; fV0Info = NULL;\r
186   if(fTracksBarrel){\r
187     fTracksBarrel->Delete(); delete fTracksBarrel;\r
188     fTracksBarrel = NULL;\r
189   }\r
190   if(fTracksSA){\r
191     fTracksSA->Delete(); delete fTracksSA;\r
192     fTracksSA = NULL;\r
193   }\r
194   if(fTracksKink){ \r
195     fTracksKink->Delete(); delete fTracksKink;\r
196     fTracksKink = NULL;\r
197   }\r
198   if(fV0List){ \r
199     fV0List->Delete(); \r
200     delete fV0List;\r
201     fV0List = NULL;\r
202   }\r
203   if(fClusters){\r
204     fClusters->Delete(); delete fClusters;\r
205     fClusters = NULL;\r
206   }\r
207   if(fContainer && !(AliAnalysisManager::GetAnalysisManager() && AliAnalysisManager::GetAnalysisManager()->IsProofMode())){\r
208     fContainer->Delete(); \r
209     delete fContainer;\r
210     fContainer = NULL;\r
211   }\r
212 }\r
213 \r
214 //____________________________________________________________________\r
215 Bool_t AliTRDinfoGen::GetRefFigure(Int_t)\r
216 {\r
217 // General graphs for PWG1/TRD train\r
218   if(!gPad){\r
219     AliWarning("Please provide a canvas to draw results.");\r
220     return kFALSE;\r
221   }\r
222   fContainer->At(kStatTrk)->Draw("bar");\r
223   return kTRUE;\r
224 }\r
225 \r
226 //____________________________________________________________________\r
227 void AliTRDinfoGen::UserCreateOutputObjects()\r
228 {       \r
229   //\r
230   // Create Output Containers (TObjectArray containing 1D histograms)\r
231   //\r
232  \r
233   fTrackInfo = new AliTRDtrackInfo();\r
234   fEventInfo = new AliTRDeventInfo();\r
235   fV0Info    = new AliTRDv0Info();\r
236   fTracksBarrel = new TObjArray(200); fTracksBarrel->SetOwner(kTRUE);\r
237   fTracksSA = new TObjArray(20); fTracksSA->SetOwner(kTRUE);\r
238   fTracksKink = new TObjArray(20); fTracksKink->SetOwner(kTRUE);\r
239   fV0List = new TObjArray(10); fV0List->SetOwner(kTRUE);\r
240   fClusters = new TObjArray(AliTRDgeometry::kNdet); fClusters->SetOwner(kTRUE);\r
241 \r
242   // define general monitor\r
243   fContainer = new TObjArray(kNclasses); fContainer->SetOwner(kTRUE);\r
244   TH1 *h=new TH1I("hStat", "Run statistics;Observable;Entries", Int_t(kNObjects), -0.5, Float_t(kNObjects)-0.5);\r
245   TAxis *ax(h->GetXaxis());\r
246   ax->SetBinLabel(Int_t(kTracksESD) + 1, "ESD");\r
247   ax->SetBinLabel(Int_t(kTracksMC) + 1, "MC");\r
248   ax->SetBinLabel(Int_t(kV0) + 1, "V0");\r
249   ax->SetBinLabel(Int_t(kTPC) + 1, "TPC");\r
250   ax->SetBinLabel(Int_t(kTRDin) + 1, "TRDin");\r
251   ax->SetBinLabel(Int_t(kTRDout) + 1, "TRDout");\r
252   ax->SetBinLabel(Int_t(kBarrel) + 1, "Barrel");\r
253   ax->SetBinLabel(Int_t(kBarrelMC) + 1, "BarrelMC");\r
254   ax->SetBinLabel(Int_t(kSA) + 1, "SA");\r
255   ax->SetBinLabel(Int_t(kSAMC) + 1, "SAMC");\r
256   ax->SetBinLabel(Int_t(kKink) + 1, "Kink");\r
257   ax->SetBinLabel(Int_t(kKinkMC) + 1, "KinkMC");\r
258   ax->SetBinLabel(Int_t(kBarrelFriend) + 1, "BFriend");\r
259   ax->SetBinLabel(Int_t(kSAFriend) + 1, "SFriend");\r
260   fContainer->AddAt(h, kStatTrk);\r
261   h=new TH1I("hEv", "Run statistics;Event Class;Entries", 4, -0.5, 3.5);\r
262   ax = h->GetXaxis();\r
263   ax->SetBinLabel(1, "Low");\r
264   ax->SetBinLabel(2, "High");\r
265   ax->SetBinLabel(3, "Cosmic");\r
266   ax->SetBinLabel(4, "Calib");\r
267   fContainer->AddAt(h, kEvType);\r
268   TH2I* h2=new TH2I("hBCtrack", "Track Statistics;Fill Bunch;TOF BC;Entries", 3500, -0.5, 3499.5, 31, -10.5, 20.5);\r
269   fContainer->AddAt(h2, kBC);\r
270   h=new TH1I("hTriggers", "Triggers statistics;;Entries", 21, -0.5, 20.5);\r
271   fContainer->AddAt(h, kTrigger);\r
272   TObjArray *chmb = new TObjArray(AliTRDgeometry::kNdet);\r
273   chmb->SetName("Chambers"); chmb->SetOwner();\r
274   fContainer->AddAt(chmb, kChmb);\r
275 \r
276   PostData(AliTRDpwg1Helper::kTracksBarrel, fTracksBarrel);\r
277   PostData(AliTRDpwg1Helper::kTracksSA,     fTracksSA);\r
278   PostData(AliTRDpwg1Helper::kTracksKink,   fTracksKink);\r
279   PostData(AliTRDpwg1Helper::kEventInfo,    fEventInfo);\r
280   PostData(AliTRDpwg1Helper::kV0List,       fV0List);\r
281   PostData(AliTRDpwg1Helper::kClusters,     fClusters);\r
282   PostData(AliTRDpwg1Helper::kMonitor,      fContainer);\r
283 }\r
284 \r
285 //____________________________________________________________________\r
286 Bool_t AliTRDinfoGen::Load(const Char_t *file, const Char_t *dir, const Char_t *name)\r
287 {\r
288 // Load data from performance file\r
289 \r
290   if(!TFile::Open(file)){\r
291     AliWarning(Form("Couldn't open file %s.", file));\r
292     return kFALSE;\r
293   }\r
294   if(dir){\r
295     if(!gFile->cd(dir)){\r
296       AliWarning(Form("Couldn't cd to %s in %s.", dir, file));\r
297       return kFALSE;\r
298     }\r
299   }\r
300   TObjArray *o(NULL);\r
301   const Char_t *tn=(name ? name : GetName());\r
302   if(!(o = (TObjArray*)gDirectory->Get(tn))){\r
303     AliWarning(Form("Missing histogram container %s.", tn));\r
304     return kFALSE;\r
305   }\r
306   fContainer = (TObjArray*)o->Clone(GetName());\r
307   gFile->Close();\r
308   return kTRUE;\r
309 }\r
310 \r
311 //____________________________________________________________________\r
312 void AliTRDinfoGen::UserExec(Option_t *){\r
313   //\r
314   // Run the Analysis\r
315   //\r
316 \r
317   fTracksBarrel->Delete();\r
318   fTracksSA->Delete();\r
319   fTracksKink->Delete();\r
320   fV0List->Delete();\r
321   fClusters->Delete();\r
322   fEventInfo->Delete("");\r
323 \r
324   fESDev = dynamic_cast<AliESDEvent*>(InputEvent());\r
325   if(!fESDev){\r
326     AliError("Failed retrieving ESD event");\r
327     return;\r
328   }\r
329   // WARNING\r
330   // This part may conflict with other detectors !!\r
331   if(!IsInitOCDB()){\r
332     AliCDBEntry* obj(NULL);\r
333     AliCDBManager* ocdb = AliCDBManager::Instance();\r
334     if(ocdb->IsDefaultStorageSet()){\r
335       AliInfo("OCDB :: initialized externally.");\r
336     } else {\r
337       AliInfo("OCDB :: initializing locally ...");\r
338       // prepare OCDB access\r
339       ocdb->SetDefaultStorage(fOCDB.Data());\r
340       ocdb->SetRun(fESDev->GetRunNumber());\r
341       // create geo manager\r
342       if(!(obj = ocdb->Get(AliCDBPath("GRP", "Geometry", "Data")))){\r
343         AliError("GEOMETRY failed initialization.");\r
344       } else {\r
345         AliGeomManager::SetGeometry((TGeoManager*)obj->GetObject());\r
346         AliGeomManager::GetNalignable("TRD");\r
347         AliGeomManager::ApplyAlignObjsFromCDB("TRD");\r
348       }\r
349       //init magnetic field\r
350       if(!TGeoGlobalMagField::Instance()->IsLocked() &&\r
351         !fESDev->InitMagneticField()){\r
352         AliError("MAGNETIC FIELD failed initialization.");\r
353       }\r
354       // set no of time bins\r
355       AliTRDtrackerV1::SetNTimeBins(AliTRDcalibDB::Instance()->GetNumberOfTimeBinsDCS());\r
356     }\r
357     AliInfo(Form("OCDB :  Loc[%s] Run[%d] TB[%d]", ocdb->GetDefaultStorage()->GetURI().Data(), ocdb->GetRun(), AliTRDtrackerV1::GetNTimeBins()));\r
358 \r
359     // load misalignment\r
360     fgGeo = new AliTRDgeometry;\r
361     fgGeo->CreateClusterMatrixArray();\r
362     MakeChambers();\r
363     // load reco param list from OCDB\r
364     AliInfo("Initializing TRD reco params ...");\r
365     fgReconstructor = new AliTRDReconstructor();\r
366     if(!(obj = ocdb->Get(AliCDBPath("TRD", "Calib", "RecoParam")))){\r
367       AliError("RECO PARAM failed initialization.");\r
368     } else {\r
369       obj->PrintMetaData();\r
370       fRecos = (TObjArray*)obj->GetObject();\r
371     }\r
372     SetInitOCDB();\r
373   }\r
374 \r
375   AliDebug(2, "+++++++++++++++++++++++++++++++++++++++");\r
376   AliDebug(2, Form("Analyse Ev[%d]", fESDev->GetEventNumberInFile()));\r
377 \r
378   // set reco param valid for this event\r
379   TH1 *h = (TH1I*)fContainer->At(kEvType);\r
380   if(!fRecos){\r
381     fgReconstructor->SetRecoParam(AliTRDrecoParam::GetLowFluxParam());\r
382     h->Fill(0);\r
383   } else {\r
384     for(Int_t ireco(0); ireco<fRecos->GetEntriesFast(); ireco++){\r
385       AliTRDrecoParam *reco((AliTRDrecoParam*)fRecos->At(ireco));\r
386       Int_t es(reco->GetEventSpecie());\r
387       if(!(es&fESDev->GetEventSpecie())) continue;\r
388       fgReconstructor->SetRecoParam(reco);\r
389       if(AliLog::GetDebugLevel("PWG1/TRD", "AliTRDinfoGen")>5) reco->Dump();\r
390       TString s;\r
391       if(es&AliRecoParam::kLowMult){ s="LowMult"; h->Fill(0);}\r
392       else if(es&AliRecoParam::kHighMult){ s="HighMult"; h->Fill(1);}\r
393       else if(es&AliRecoParam::kCosmic){ s="Cosmic"; h->Fill(2);}\r
394       else if(es&AliRecoParam::kCalib){ s="Calib"; h->Fill(3);}\r
395       else s="Unknown";\r
396       AliDebug(2, Form("  Using reco param \"%s\" for event %d.", s.Data(), fESDev->GetEventNumberInFile()));\r
397       break;\r
398     }\r
399   }\r
400 \r
401   // link MC if available\r
402   fMCev = MCEvent();\r
403   \r
404   // trigger monitor\r
405   h = (TH1I*)fContainer->At(kTrigger);\r
406   TAxis *ax(h->GetXaxis());\r
407   TObjArray *evTriggers = fESDev->GetFiredTriggerClasses().Tokenize(" ");\r
408   for(Int_t iet(evTriggers->GetEntriesFast()); iet--;){\r
409     Int_t ix(1);\r
410     for(; ix<=ax->GetNbins(); ix++){\r
411       if(!Int_t(h->GetBinContent(ix))){\r
412         ax->SetBinLabel(ix, (*evTriggers)[iet]->GetName());\r
413         break;\r
414       }\r
415       if(strcmp((*evTriggers)[iet]->GetName(), ax->GetBinLabel(ix))==0) break;\r
416     }\r
417     h->AddBinContent(ix);\r
418   }\r
419 \r
420   // event selection based on vertex cuts and trigger\r
421   if(UseLocalEvSelection() && !fEventCut->IsSelected(fESDev, IsCollision())) return;\r
422 \r
423   if(!fESDfriend){\r
424     AliError("Failed retrieving ESD friend event");\r
425     return;\r
426   }\r
427   if(HasMCdata() && !fMCev){\r
428     AliError("Failed retrieving MC event");\r
429     return;\r
430   }\r
431 \r
432   new(fEventInfo)AliTRDeventInfo(fESDev->GetHeader(), const_cast<AliESDRun *>(fESDev->GetESDRun()));\r
433   // Determine centrality\r
434   // Author: Ionut Arsene <I.C.Arsene@gsi.de>\r
435   AliDebug(2, Form("  Beam Type: %s", fESDev->GetESDRun()->GetBeamType()));\r
436   TString beamtype = fESDev->GetESDRun()->GetBeamType();\r
437   if(beamtype.Contains("Pb-Pb") || beamtype.Contains("A-A")){\r
438     const AliMultiplicity *mult = fESDev->GetMultiplicity();\r
439     fEventInfo->SetMultiplicity(mult?mult->GetNumberOfTracklets():0);\r
440     const AliCentrality *cent = fESDev->GetCentrality();\r
441     // centrality for different options V0 = "V0M", ITS = "TKL" etc\r
442     fEventInfo->SetCentrality(cent?cent->GetCentralityPercentile("TKL"):-1.);\r
443   } else {\r
444     fEventInfo->SetMultiplicity(0);\r
445     fEventInfo->SetCentrality(-1.);\r
446   }\r
447   UShort_t evBC(fESDev->GetBunchCrossNumber());\r
448 \r
449   Bool_t *trackMap(NULL);\r
450   AliStack * mStack(NULL);\r
451   Int_t nTracksMC = HasMCdata() ? fMCev->GetNumberOfTracks() : 0, nTracksESD = fESDev->GetNumberOfTracks();\r
452   if(HasMCdata()){\r
453     mStack = fMCev->Stack();\r
454     if(!mStack){\r
455       AliError("Failed retrieving MC Stack");\r
456       return;\r
457     }\r
458     trackMap = new Bool_t[nTracksMC];\r
459     memset(trackMap, 0, sizeof(Bool_t) * nTracksMC);\r
460   }\r
461   \r
462   Double32_t dedx[100]; Int_t nSlices(0);\r
463   Int_t nTRDout(0), nTRDin(0), nTPC(0)\r
464        ,nclsTrklt\r
465        ,nBarrel(0), nSA(0), nKink(0)\r
466        ,nBarrelFriend(0), nSAFriend(0)\r
467        ,nBarrelMC(0), nSAMC(0), nKinkMC(0);\r
468   AliESDtrack *esdTrack = NULL;\r
469   AliESDfriendTrack *esdFriendTrack = NULL;\r
470   TObject *calObject = NULL;\r
471   AliTRDtrackV1 *track = NULL;\r
472   AliTRDseedV1 *tracklet = NULL;\r
473   AliTRDcluster *cl = NULL;\r
474 \r
475 \r
476   // LOOP 0 - over ESD v0s\r
477   Float_t bField(fESDev->GetMagneticField());\r
478   AliESDv0 *v0(NULL);\r
479   Int_t v0pid[AliPID::kSPECIES];\r
480   for(Int_t iv0(0); iv0<fESDev->GetNumberOfV0s(); iv0++){\r
481     if(!(v0 = fESDev->GetV0(iv0))) continue;\r
482     // register v0\r
483     if(fV0Cut) new(fV0Info) AliTRDv0Info(*fV0Cut);\r
484     else new(fV0Info) AliTRDv0Info();\r
485     fV0Info->SetMagField(bField);\r
486     fV0Info->SetV0tracks(fESDev->GetTrack(v0->GetPindex()), fESDev->GetTrack(v0->GetNindex()));\r
487     fV0Info->SetV0Info(v0);\r
488     fV0List->Add(new AliTRDv0Info(*fV0Info));//  kFOUND=kFALSE;\r
489   }\r
490 \r
491 \r
492   // LOOP 1 - over ESD tracks\r
493   AliTRDv0Info *v0info=NULL;\r
494   for(Int_t itrk = 0; itrk < nTracksESD; itrk++){\r
495     new(fTrackInfo) AliTRDtrackInfo();\r
496     esdTrack = fESDev->GetTrack(itrk);\r
497     AliDebug(3, Form("\n%3d ITS[%d] TPC[%d] TRD[%d] TOF-BC[%d]\n", itrk, esdTrack->GetNcls(0), esdTrack->GetNcls(1), esdTrack->GetNcls(2), esdTrack->GetTOFBunchCrossing()));\r
498     if(esdTrack->GetStatus()&AliESDtrack::kTPCout) nTPC++;\r
499     if(esdTrack->GetStatus()&AliESDtrack::kTRDout) nTRDout++;\r
500     if(esdTrack->GetStatus()&AliESDtrack::kTRDin) nTRDin++;\r
501     // look at external track param\r
502     const AliExternalTrackParam *op = esdTrack->GetOuterParam();\r
503     Double_t xyz[3];\r
504     if(op){\r
505       op->GetXYZ(xyz);\r
506       op->Global2LocalPosition(xyz, op->GetAlpha());\r
507       AliDebug(3, Form("op @ X[%7.3f]\n", xyz[0]));\r
508     }\r
509 \r
510     // read MC info\r
511     Int_t fPdg = -1;\r
512     Int_t label = -1; UInt_t alab=UINT_MAX;\r
513     if(HasMCdata()){\r
514       label = esdTrack->GetLabel(); \r
515       alab = TMath::Abs(label);\r
516       // register the track\r
517       if(alab < UInt_t(nTracksMC)){ \r
518         trackMap[alab] = kTRUE; \r
519       } else { \r
520         AliError(Form("MC label[%d] outside scope for Ev[%d] Trk[%d].", label, (Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), itrk));\r
521         continue; \r
522       }\r
523       AliMCParticle *mcParticle = NULL; \r
524       if(!(mcParticle = (AliMCParticle*) fMCev->GetTrack(alab))){\r
525         AliError(Form("MC particle label[%d] missing for Ev[%d] Trk[%d].", label, (Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), itrk));\r
526         continue;\r
527       }\r
528       fPdg = mcParticle->Particle()->GetPdgCode();\r
529       Int_t nRefs = mcParticle->GetNumberOfTrackReferences();\r
530       Int_t iref = 0; AliTrackReference *ref = NULL; \r
531       while(iref<nRefs){\r
532         ref = mcParticle->GetTrackReference(iref);\r
533         if(ref->LocalX() > fgkTPC) break;\r
534         iref++;\r
535       }\r
536 \r
537       fTrackInfo->SetMC();\r
538       fTrackInfo->SetPDG(fPdg);\r
539       fTrackInfo->SetPrimary(mcParticle->Particle()->IsPrimary());\r
540       fTrackInfo->SetLabel(label);\r
541       Int_t jref = iref;//, kref = 0;\r
542       while(jref<nRefs){\r
543         ref = mcParticle->GetTrackReference(jref);\r
544         if(ref->LocalX() > fgkTRD) break;\r
545         AliDebug(4, Form("  trackRef[%2d (%2d)] @ %7.3f OK", jref-iref, jref, ref->LocalX()));\r
546         fTrackInfo->AddTrackRef(ref);\r
547         jref++;\r
548       }\r
549       AliDebug(3, Form("NtrackRefs[%d(%d)]", fTrackInfo->GetNTrackRefs(), nRefs));\r
550     }\r
551 \r
552     // copy some relevant info to TRD track info\r
553     fTrackInfo->SetStatus(esdTrack->GetStatus());\r
554     fTrackInfo->SetTrackId(esdTrack->GetID());\r
555     Double_t p[AliPID::kSPECIES]; esdTrack->GetTRDpid(p);\r
556     fTrackInfo->SetESDpid(p);\r
557     fTrackInfo->SetESDpidQuality(esdTrack->GetTRDntrackletsPID());\r
558     if(!nSlices) nSlices = esdTrack->GetNumberOfTRDslices();\r
559     memset(dedx, 0, 100*sizeof(Double32_t));\r
560     Int_t in(0);\r
561     for(Int_t il=0; il<AliTRDgeometry::kNlayer; il++)\r
562       for(Int_t is=0; is<nSlices; is++) \r
563         dedx[in++]=esdTrack->GetTRDslice(il, is);\r
564     for(Int_t il=0; il<AliTRDgeometry::kNlayer; il++) dedx[in++]=esdTrack->GetTRDmomentum(il);\r
565     fTrackInfo->SetSlices(in, dedx);\r
566     fTrackInfo->SetNumberOfClustersRefit(esdTrack->GetNcls(2));\r
567     // some other Informations which we may wish to store in order to find problematic cases\r
568     fTrackInfo->SetKinkIndex(esdTrack->GetKinkIndex(0));\r
569     fTrackInfo->SetTPCncls(static_cast<UShort_t>(esdTrack->GetNcls(1)));\r
570     fTrackInfo->SetTOFbc(esdTrack->GetTOFBunchCrossing()==AliVTrack::kTOFBCNA?0:esdTrack->GetTOFBunchCrossing());\r
571     nclsTrklt = 0;\r
572   \r
573     // set V0pid info\r
574     for(Int_t iv(0); iv<fV0List->GetEntriesFast(); iv++){\r
575       if(!(v0info = (AliTRDv0Info*)fV0List->At(iv))) continue;\r
576       if(!v0info->GetV0Daughter(1) && !v0info->GetV0Daughter(-1)) continue;\r
577       if(!v0info->HasTrack(fTrackInfo)) continue;\r
578       memset(v0pid, 0, AliPID::kSPECIES*sizeof(Int_t));\r
579       fTrackInfo->SetV0();\r
580       for(Int_t is=AliPID::kSPECIES; is--;){v0pid[is] = v0info->GetPID(is, fTrackInfo);}\r
581       fTrackInfo->SetV0pid(v0pid);\r
582       fTrackInfo->SetV0();\r
583       //const AliTRDtrackInfo::AliESDinfo *ei = fTrackInfo->GetESDinfo();\r
584       break;\r
585     }\r
586 \r
587     // read track REC info\r
588     esdFriendTrack = (fESDfriend->GetNumberOfTracks() > itrk) ? fESDfriend->GetTrack(itrk): NULL;\r
589     if(esdFriendTrack){\r
590       fTrackInfo->SetTPCoutParam(esdFriendTrack->GetTPCOut());\r
591       Int_t icalib = 0;\r
592       while((calObject = esdFriendTrack->GetCalibObject(icalib++))){\r
593         if(strcmp(calObject->IsA()->GetName(),"AliTRDtrackV1") != 0) continue; // Look for the TRDtrack\r
594         if(!(track = dynamic_cast<AliTRDtrackV1*>(calObject))) break;\r
595         AliDebug(4, Form("TRD track OK"));\r
596         // Set the clusters to unused\r
597         for(Int_t ipl = 0; ipl < AliTRDgeometry::kNlayer; ipl++){\r
598           if(!(tracklet = track->GetTracklet(ipl))) continue;\r
599           tracklet->ResetClusterIter();\r
600           while((cl = tracklet->NextCluster())) cl->Use(0);\r
601         }\r
602         fTrackInfo->SetTrack(track);\r
603         ((TH2I*)fContainer->At(kBC))->Fill(evBC, esdTrack->GetTOFBunchCrossing());\r
604         break;\r
605       }\r
606       AliDebug(3, Form("Ntracklets[%d]\n", fTrackInfo->GetNTracklets()));\r
607     } else AliDebug(3, "No ESD friends");\r
608     if(op) fTrackInfo->SetOuterParam(op);\r
609 \r
610     if(DebugLevel() >= 1){\r
611       AliTRDtrackInfo info(*fTrackInfo);\r
612       (*DebugStream()) << "trackInfo"\r
613       << "TrackInfo.=" << &info\r
614       << "\n";\r
615       info.Delete("");\r
616     }\r
617 \r
618     ULong_t status(esdTrack->GetStatus());\r
619     if((status&AliESDtrack::kTPCout)){\r
620       if(!esdTrack->GetKinkIndex(0)){ // Barrel  Track Selection\r
621         Bool_t selected(kTRUE);\r
622         if(UseLocalTrkSelection()){\r
623           if(esdTrack->Pt() < fgkPt){ \r
624             AliDebug(3, Form("Reject Trk[%3d] Ev[%4d] Pt[%5.2f]", itrk, fESDev->GetEventNumberInFile(), esdTrack->Pt()));\r
625             selected = kFALSE;\r
626           }\r
627           if(selected && TMath::Abs(esdTrack->Eta()) > fgkEta){\r
628             AliDebug(3, Form("Reject Trk[%3d] Ev[%4d] Eta[%5.2f]", itrk, fESDev->GetEventNumberInFile(), TMath::Abs(esdTrack->Eta())));\r
629             selected = kFALSE;\r
630           }\r
631           if(selected && esdTrack->GetTPCNcls() < fgkNclTPC){ \r
632             AliDebug(3, Form("Reject Trk[%3d] Ev[%4d] NclTPC[%d]", itrk, fESDev->GetEventNumberInFile(), esdTrack->GetTPCNcls()));\r
633             selected = kFALSE;\r
634           }\r
635           Float_t par[2], cov[3];\r
636           esdTrack->GetImpactParameters(par, cov);\r
637           if(IsCollision()){ // cuts on DCA\r
638             if(selected && TMath::Abs(par[0]) > fgkTrkDCAxy){ \r
639               AliDebug(3, Form("Reject Trk[%3d] Ev[%4d] DCAxy[%f]", itrk, fESDev->GetEventNumberInFile(), TMath::Abs(par[0])));\r
640               selected = kFALSE;\r
641             }\r
642             if(selected && TMath::Abs(par[1]) > fgkTrkDCAz){ \r
643               AliDebug(3, Form("Reject Trk[%3d] Ev[%4d] DCAz[%f]", itrk, fESDev->GetEventNumberInFile(), TMath::Abs(par[1])));\r
644               selected = kFALSE;\r
645             }\r
646           } else if(selected && fMCev && !fMCev->IsPhysicalPrimary(alab)){;\r
647             AliDebug(3, Form("Reject Trk[%3d] Ev[%4d] Primary", itrk, fESDev->GetEventNumberInFile()));\r
648             selected = kFALSE;\r
649           }\r
650         }\r
651         if(fTrackCut && !fTrackCut->IsSelected(esdTrack)) selected = kFALSE;\r
652         if(selected){ \r
653           fTracksBarrel->Add(new AliTRDtrackInfo(*fTrackInfo));\r
654           nBarrel++;\r
655           if(fTrackInfo->GetTrack()) \r
656             nBarrelFriend++;\r
657         }\r
658       } else {\r
659         fTracksKink->Add(new AliTRDtrackInfo(*fTrackInfo));\r
660         nKink++;\r
661       }\r
662     } else if((status&AliESDtrack::kTRDout) && !(status&AliESDtrack::kTRDin)){ \r
663       fTracksSA->Add(new AliTRDtrackInfo(*fTrackInfo));\r
664       nSA++;\r
665       if(fTrackInfo->GetTrack()) \r
666         nSAFriend++;\r
667     }\r
668     fTrackInfo->Delete("");\r
669   }\r
670   // read clusters REC info\r
671   TTree * treeR = (dynamic_cast<AliESDInputHandlerRP*>(fInputHandler))->GetTreeR("TRD");\r
672   if(treeR) {\r
673     TObjArray *recPoints(NULL);\r
674     if((treeR->GetBranch("TRDcluster"))){\r
675       treeR->SetBranchAddress("TRDcluster", &recPoints);\r
676       for(Int_t idet(0); idet<treeR->GetEntries(); idet++){\r
677         treeR->GetEntry(idet);\r
678         if(!recPoints->GetEntries()){\r
679           AliDebug(1, Form("Missing entry %d from TreeR", idet));\r
680           continue;\r
681         }\r
682         AliTRDcluster *c = (AliTRDcluster*)(*recPoints)[0];\r
683         if(!c){\r
684           AliDebug(1, Form("Missing first cluster in entry %d from TreeR", idet));\r
685           continue;\r
686         }\r
687         fClusters->AddAt(recPoints->Clone(Form("%03d", c->GetDetector())), c->GetDetector());\r
688       }\r
689     } else AliDebug(3, "No TRDcluster branch");\r
690   } else AliDebug(3, "No RecPoints");\r
691 \r
692 \r
693 \r
694   // LOOP 2 - over MC tracks which are passing TRD where the track is not reconstructed\r
695   if(HasMCdata()){\r
696     AliDebug(10, "Output of the MC track map:");\r
697     for(Int_t itk = 0; itk < nTracksMC;  itk++) AliDebug(10, Form("trackMap[%d] = %s", itk, trackMap[itk] == kTRUE ? "TRUE" : "kFALSE"));\r
698   \r
699     for(Int_t itk = 0; itk < nTracksMC; itk++){\r
700       if(trackMap[itk]) continue;\r
701       AliMCParticle *mcParticle =  (AliMCParticle*) fMCev->GetTrack(TMath::Abs(itk));\r
702       Int_t fPdg = mcParticle->Particle()->GetPdgCode();\r
703       Int_t nRefs = mcParticle->GetNumberOfTrackReferences();\r
704       Int_t iref = 0; AliTrackReference *ref = NULL; \r
705       Int_t nRefsTRD = 0;\r
706       new(fTrackInfo) AliTRDtrackInfo();\r
707       fTrackInfo->SetMC();\r
708       fTrackInfo->SetPDG(fPdg);\r
709       while(iref<nRefs){ // count TRD TR\r
710         Bool_t kIN(kFALSE);\r
711         ref = mcParticle->GetTrackReference(iref);\r
712         if(ref->LocalX() > fgkTPC && ref->LocalX() < fgkTRD){\r
713           fTrackInfo->AddTrackRef(ref);\r
714           nRefsTRD++;kIN=kTRUE;\r
715         }\r
716         AliDebug(4, Form("  trackRef[%2d] @ x[%7.3f] %s", iref, ref->LocalX(), kIN?"IN":"OUT"));\r
717         iref++;\r
718       }\r
719       if(!nRefsTRD){\r
720         // In this stage we at least require 1 hit inside TRD. What will be done with this tracks is a task for the \r
721         // analysis job\r
722         fTrackInfo->Delete("");\r
723         continue;\r
724       }\r
725       fTrackInfo->SetPrimary(mcParticle->Particle()->IsPrimary());\r
726       fTrackInfo->SetLabel(itk);\r
727       if(DebugLevel() >= 1){\r
728         AliTRDtrackInfo info(*fTrackInfo);\r
729         (*DebugStream()) << "trackInfo"\r
730         << "TrackInfo.=" << &info\r
731         << "\n";\r
732         info.Delete("");\r
733       }\r
734       AliDebug(3, Form("Add MC track @ label[%d] nTRDrefs[%d].", itk, nRefsTRD));\r
735       // check where the track starts\r
736       ref = mcParticle->GetTrackReference(0);\r
737       if(ref->LocalX() < fgkITS){ \r
738         fTracksBarrel->Add(new AliTRDtrackInfo(*fTrackInfo));\r
739         nBarrelMC++;\r
740       } else if(ref->LocalX() < fgkTPC) {\r
741         fTracksKink->Add(new AliTRDtrackInfo(*fTrackInfo));\r
742         nKinkMC++;\r
743       } else if(nRefsTRD>6){\r
744         fTracksSA->Add(new AliTRDtrackInfo(*fTrackInfo));\r
745         nSAMC++;\r
746       }\r
747       fTrackInfo->Delete("");\r
748     }\r
749     delete[] trackMap;\r
750   }\r
751   AliDebug(1, Form(\r
752     "\nEv[%3d] Tracks: ESD[%d] MC[%d] V0[%d]\n"\r
753     "        TPCout[%d] TRDin[%d] TRDout[%d]\n"\r
754     "        Barrel[%3d+%3d=%3d] SA[%2d+%2d=%2d] Kink[%2d+%2d=%2d]"\r
755     ,(Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), nTracksESD, nTracksMC, fV0List->GetEntries()\r
756     , nTPC, nTRDin, nTRDout\r
757     ,nBarrel, nBarrelMC, fTracksBarrel->GetEntries()\r
758     ,nSA, nSAMC, fTracksSA->GetEntries()\r
759     ,nKink, nKinkMC, fTracksKink->GetEntries()\r
760   ));\r
761   // save track statistics\r
762   h = (TH1I*)fContainer->At(kStatTrk);\r
763   h->Fill(Float_t(kTracksESD), nTracksESD);\r
764   h->Fill(Float_t(kTracksMC), nTracksMC);\r
765   h->Fill(Float_t(kV0), fV0List->GetEntries());\r
766   h->Fill(Float_t(kTPC), nTPC);\r
767   h->Fill(Float_t(kTRDin), nTRDin);\r
768   h->Fill(Float_t(kTRDout), nTRDout);\r
769   h->Fill(Float_t(kBarrel), nBarrel);\r
770   h->Fill(Float_t(kBarrelMC), nBarrelMC);\r
771   h->Fill(Float_t(kSA), nSA);\r
772   h->Fill(Float_t(kSAMC), nSAMC);\r
773   h->Fill(Float_t(kKink), nKink);\r
774   h->Fill(Float_t(kKinkMC), nKinkMC);\r
775   h->Fill(Float_t(kBarrelFriend), nBarrelFriend);\r
776   h->Fill(Float_t(kSAFriend), nSAFriend);\r
777 }\r
778 \r
779 \r
780 //____________________________________________________________________\r
781 void AliTRDinfoGen::MakeChambers()\r
782 {\r
783 // Build chamber position and status\r
784   if(!fContainer){\r
785     AliError("Missing container");\r
786     return;\r
787   }\r
788   AliTRDcalibDB *calib(AliTRDcalibDB::Instance());\r
789   if(!calib){\r
790     AliError("No access to calibration data");\r
791     return;\r
792   }\r
793 \r
794   TObjArray *chmb((TObjArray*)fContainer->At(kChmb));\r
795   Double_t alpha(0.), cs(-2.), sn(0.);\r
796   TVectorF pos(5);\r
797   for(Int_t isec(0); isec<AliTRDgeometry::kNsector; isec++){\r
798     alpha = (0.5+isec)*AliTRDgeometry::GetAlpha();\r
799     cs    = TMath::Cos(alpha);\r
800     sn    = TMath::Sin(alpha);\r
801 \r
802     for(Int_t istk(0); istk<AliTRDgeometry::kNstack; istk++){\r
803       for(Int_t ilyr(0); ilyr<AliTRDgeometry::kNlayer; ilyr++){\r
804         Int_t idet(AliTRDgeometry::GetDetector(ilyr, istk, isec));\r
805         TGeoHMatrix *matrix(fgGeo->GetClusterMatrix(idet));\r
806         if(!matrix){\r
807           AliDebug(2, Form("Missing matrix for %03d [%02d_%d_%d]", idet, isec, istk, ilyr));\r
808           continue;\r
809         }\r
810         AliDebug(2, Form("Read info for %03d [%02d_%d_%d]", idet, isec, istk, ilyr));\r
811         AliTRDpadPlane *pp(fgGeo->GetPadPlane(ilyr, istk));\r
812         Double_t zm(0.5 * (pp->GetRow0() + pp->GetRowEnd())),\r
813                  loc0[] = {AliTRDgeometry::AnodePos(), pp->GetCol0(), zm-pp->GetRow0()},\r
814                  loc1[] = {AliTRDgeometry::AnodePos(), pp->GetColEnd(), zm-pp->GetRowEnd()},\r
815                  glb[3];\r
816         matrix->LocalToMaster(loc0, glb);\r
817         Float_t phi = TMath::ATan2(glb[0]*sn + glb[1]*cs, glb[0]*cs - glb[1]*sn),\r
818                 tgl = glb[2]/glb[0]/TMath::Sqrt(1.+glb[1]*glb[1]/glb[0]/glb[0]),\r
819                 eta = -TMath::Log(TMath::Tan(0.5 *  (0.5*TMath::Pi() - TMath::ATan(tgl))));\r
820         pos[0] = eta; pos[1] = phi;\r
821         matrix->LocalToMaster(loc1, glb);\r
822         phi = TMath::ATan2(glb[0]*sn + glb[1]*cs, glb[0]*cs - glb[1]*sn);\r
823         tgl = glb[2]/glb[0]/TMath::Sqrt(1.+glb[1]*glb[1]/glb[0]/glb[0]);\r
824         eta = -TMath::Log(TMath::Tan(0.5 *  (0.5*TMath::Pi() - TMath::ATan(tgl))));\r
825         pos[2] = eta; pos[3] = phi;\r
826         pos[4] = 0.;\r
827         if(calib->IsChamberGood(idet)){\r
828           if(calib->IsHalfChamberNoData(idet, 0)) pos[4] += 2.;\r
829           if(calib->IsHalfChamberNoData(idet, 1)) pos[4] += 3.;\r
830         } else pos[4] = 1.;\r
831         chmb->AddAt(new TVectorF(pos), idet);\r
832       }\r
833     }\r
834   }\r
835 }\r
836 \r
837 //____________________________________________________________________\r
838 void AliTRDinfoGen::MakeSummary()\r
839 {\r
840 // Build summary plots\r
841   if(!fContainer){\r
842     AliError("Missing results");\r
843     return;\r
844   }\r
845   TH1 *h1(NULL); TVirtualPad *p(NULL); TCanvas *cOut(NULL);\r
846 \r
847   const Int_t nx(2048), ny(750);\r
848   cOut = new TCanvas(GetName(), "Run Statistics", nx, ny);\r
849   cOut->Divide(3,1, 1.e-5, 1.e-5);\r
850   //=========\r
851   p=cOut->cd(1);p->SetRightMargin(0.025);p->SetTopMargin(0.01);p->SetLogy();\r
852   h1 = (TH1*)fContainer->At(kStatTrk);\r
853   h1->SetBarOffset(0.06); h1->SetBarWidth(0.88); h1->SetFillColor(3);\r
854   h1->Draw("bar1");\r
855   //=========\r
856   p=cOut->cd(2);p->SetRightMargin(0.025);p->SetTopMargin(0.01);p->SetLogy();\r
857   h1 = (TH1*)fContainer->At(kEvType);\r
858   h1->SetBarOffset(0.04); h1->SetBarWidth(0.92);h1->SetFillColor(6);\r
859   h1->Draw("bar1");\r
860   //=========\r
861   p=cOut->cd(3);p->SetRightMargin(0.025);p->SetTopMargin(0.01);p->SetLogz();\r
862   TH2 *h2((TH2*)fContainer->At(kBC));\r
863   h1 = h2->ProjectionX("_px");\r
864   Int_t n(0); Int_t bins[3500];\r
865   for(Int_t ib(1); ib<=3500; ib++){\r
866     if(h1->GetBinContent(ib) < 1) continue;\r
867     bins[n++] = ib;\r
868   }\r
869   delete h1;\r
870 \r
871   TAxis *ay(h2->GetYaxis());\r
872   TH2 *hs = new TH2I("hBC_Summary", Form("%s;%s;%s", h2->GetTitle(), h2->GetXaxis()->GetTitle(), h2->GetYaxis()->GetTitle()),\r
873                      n, -0.5, n-0.5, ay->GetNbins(), ay->GetXmin(), ay->GetXmax());\r
874   hs->SetLineColor(kBlack);hs->SetLineWidth(1);\r
875   hs->SetMarkerColor(kRed);\r
876   TAxis *ax(hs->GetXaxis());\r
877   for(Int_t ib(0); ib<n; ib++){\r
878     ax->SetBinLabel(ib+1, Form("%d", bins[ib]));\r
879     for(Int_t iy(1); iy<=ay->GetNbins(); iy++){\r
880       hs->SetBinContent(ib+1, iy, h2->GetBinContent(bins[ib], iy));\r
881     }\r
882   }\r
883   hs->Draw("textbox");\r
884   cOut->SaveAs(Form("%s.gif", cOut->GetName()));\r
885 }\r
886 \r
887 //____________________________________________________________________\r
888 void AliTRDinfoGen::SetLocalEvSelection(const AliTRDeventCuts &ec)\r
889 {\r
890 // Set event cuts from outside\r
891   if(!fEventCut) fEventCut = new AliTRDeventCuts(ec);\r
892   else new(fEventCut) AliTRDeventCuts(ec);\r
893   fEventCut->Print();\r
894 }\r
895 \r
896 //____________________________________________________________________\r
897 void AliTRDinfoGen::SetLocalV0Selection(const AliTRDv0Info &v0)\r
898 {\r
899 // Set V0 cuts from outside\r
900 \r
901   if(!fV0Cut) fV0Cut = new AliTRDv0Info(v0);\r
902   else new(fV0Cut) AliTRDv0Info(v0);\r
903   fV0Cut->Print();\r
904 }\r
905 \r
906 //____________________________________________________________________\r
907 TTreeSRedirector* AliTRDinfoGen::DebugStream()\r
908 {\r
909 // Manage debug stream for task\r
910   if(!fDebugStream){\r
911     TDirectory *savedir = gDirectory;\r
912     fDebugStream = new TTreeSRedirector("TRD.DebugInfoGen.root");\r
913     savedir->cd();\r
914   }\r
915   return fDebugStream;\r
916 }\r
917 \r
918 //____________________________________________________________________\r
919 void AliTRDinfoGen::Terminate(Option_t* /*option*/)\r
920 {\r
921 // Process run information\r
922   AliInfo("");\r
923   if(!(fContainer = dynamic_cast<TObjArray *>(GetOutputData(AliTRDpwg1Helper::kMonitor)))) return;\r
924   AliInfo(Form("fContainer(%p)", (void*)fContainer));\r
925   if(UseLocalEvSelection()){\r
926     TH1 *h1 = (TH1*)fContainer->At(kTrigger); TAxis *ax(h1->GetXaxis());\r
927     AliInfo(Form("h1(%p)", (void*)h1));\r
928     for(Int_t ix(1); ix<=ax->GetNbins(); ix++){\r
929       if(fEventCut->CheckTrigger(ax->GetBinLabel(ix))) ax->SetBinLabel(ix, Form("#color[2]{%s}", ax->GetBinLabel(ix)));\r
930     }\r
931   }\r
932 }\r
933 \r