]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/TRD/AliTRDinfoGen.cxx
- call UserExec instead of Exec of base class AliTRDreco
[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 <TH1F.h>
38 #include <TFile.h>
39 #include <TTree.h>
40 #include <TROOT.h>
41 #include <TChain.h>
42 #include <TParticle.h>
43
44 #include "AliLog.h"
45 #include "AliAnalysisManager.h"
46 #include "AliESDEvent.h"
47 #include "AliMCEvent.h"
48 #include "AliESDInputHandler.h"
49 #include "AliMCEventHandler.h"
50
51 #include "AliESDfriend.h"
52 #include "AliESDfriendTrack.h"
53 #include "AliESDHeader.h"
54 #include "AliESDtrack.h"
55 #include "AliESDtrackCuts.h"
56 #include "AliMCParticle.h"
57 #include "AliPID.h"
58 #include "AliStack.h"
59 #include "AliTRDtrackV1.h"
60 #include "AliTrackReference.h"
61 #include "AliTRDgeometry.h"
62 #include "AliTRDcluster.h"
63 #include "AliTRDseedV1.h"
64 #include "TTreeStream.h"
65
66 #include <cstdio>
67 #include <climits>
68 #include <cstring>
69 #include <iostream>
70
71 #include "AliTRDinfoGen.h"
72 #include "info/AliTRDtrackInfo.h"
73 #include "info/AliTRDeventInfo.h"
74 #include "info/AliTRDv0Info.h"
75 #include "info/AliTRDeventCuts.h"
76
77 ClassImp(AliTRDinfoGen)
78
79 const Float_t AliTRDinfoGen::fgkTPC = 290.;
80 const Float_t AliTRDinfoGen::fgkTOF = 365.;
81
82 const Float_t AliTRDinfoGen::fgkEvVertexZ = 15.;
83 const Int_t   AliTRDinfoGen::fgkEvVertexN = 1;
84
85 const Float_t AliTRDinfoGen::fgkTrkDCAxy  = 40.;
86 const Float_t AliTRDinfoGen::fgkTrkDCAz   = 15.;
87 const Int_t   AliTRDinfoGen::fgkNclTPC    = 10;
88 const Float_t AliTRDinfoGen::fgkPt        = 0.2;
89 const Float_t AliTRDinfoGen::fgkEta       = 0.9;
90
91 //____________________________________________________________________
92 AliTRDinfoGen::AliTRDinfoGen()
93   :AliTRDrecoTask()
94   ,fEvTrigger("CINT1B-ABCE-NOPF-ALL")
95   ,fESDev(NULL)
96   ,fMCev(NULL)
97   ,fTrackInfo(NULL)
98   ,fEventInfo(NULL)
99   ,fV0container(NULL)
100   ,fV0Info(NULL)
101   ,fEventCut(NULL)
102   ,fTrackCut(NULL)
103 {
104   //
105   // Default constructor
106   //
107   SetNameTitle("infoGen", "MC-REC TRD-track list generator");
108 }
109
110 //____________________________________________________________________
111 AliTRDinfoGen::AliTRDinfoGen(char* name)
112   :AliTRDrecoTask(name)
113   ,fEvTrigger("")
114   ,fESDev(NULL)
115   ,fMCev(NULL)
116   ,fTrackInfo(NULL)
117   ,fEventInfo(NULL)
118   ,fV0container(NULL)
119   ,fV0Info(NULL)
120   ,fEventCut(NULL)
121   ,fTrackCut(NULL)
122 {
123   //
124   // Default constructor
125   //
126     DefineOutput(2, AliTRDeventInfo::Class());   // -> TRDeventInfo
127     DefineOutput(3, TObjArray::Class());         // -> TObjArray
128 }
129
130 //____________________________________________________________________
131 AliTRDinfoGen::~AliTRDinfoGen()
132 {
133 // Destructor
134
135   if(fTrackCut) delete fTrackCut;
136   if(fEventCut) delete fEventCut;
137   if(fTrackInfo) delete fTrackInfo; fTrackInfo = NULL;
138   if(fEventInfo) delete fEventInfo; fEventInfo = NULL;
139   if(fV0Info) delete fV0Info; fV0Info = NULL;
140   if(fContainer){ 
141     fContainer->Delete(); delete fContainer;
142     fContainer = NULL;
143   }
144   if(fV0container){ 
145     fV0container->Delete(); delete fV0container;
146     fV0container = NULL;
147   }
148 }
149
150 //____________________________________________________________________
151 void AliTRDinfoGen::UserCreateOutputObjects()
152 {       
153   //
154   // Create Output Containers (TObjectArray containing 1D histograms)
155   //
156  
157   fTrackInfo = new AliTRDtrackInfo();
158   fEventInfo = new AliTRDeventInfo();
159   fV0Info    = new AliTRDv0Info();
160   fContainer = new TObjArray(1000);
161   fContainer->SetOwner(kTRUE);
162   fV0container = new TObjArray(50);
163   fV0container->SetOwner(kTRUE);
164 }
165
166 //____________________________________________________________________
167 void AliTRDinfoGen::UserExec(Option_t *){
168   //
169   // Run the Analysis
170   //
171   printf("AliTRDinfoGen::UserExec : ESD[%p] MC[%p]\n", (void*)fESDev, (void*)fMCev);
172
173   fESDev = dynamic_cast<AliESDEvent*>(InputEvent());
174   fMCev = MCEvent();
175   AliInfo(Form("ESD[%p] MC[%p]", (void*)fESDev, (void*)fMCev));
176
177   if(!fESDev){
178     AliError("Failed retrieving ESD event");
179     return;
180   }
181
182   // event selection : trigger cut
183   if(UseLocalEvSelection() &&
184      fEvTrigger.Data()[0]!='\0'){ 
185     Bool_t kTRIGGERED(kFALSE);
186     const TObjArray *trig = fEvTrigger.Tokenize(" ");
187     for(Int_t itrig=trig->GetEntriesFast(); itrig--;){
188       const Char_t *trigClass(((TObjString*)(*trig)[itrig])->GetName());
189       if(fESDev->IsTriggerClassFired(trigClass)) {
190         AliDebug(2, Form("Ev[%4d] Trigger[%s]", fESDev->GetEventNumberInFile(), trigClass));
191         kTRIGGERED = kTRUE;
192         break; 
193       }
194     }
195     if(!kTRIGGERED) return;
196
197     // select only physical events
198     if(fESDev->GetEventType() != 7) return;
199   }
200
201   // if the required trigger is a collision trigger then apply event vertex cut
202   if(UseLocalEvSelection() && IsCollision()){
203     const AliESDVertex *vertex = fESDev->GetPrimaryVertex();
204     if(TMath::Abs(vertex->GetZv())<1.e-10 || 
205        TMath::Abs(vertex->GetZv())>fgkEvVertexZ || 
206        vertex->GetNContributors()<fgkEvVertexN) {
207       //PostData(0, fOutStorage);
208       return;
209     }
210   }
211
212   if(fEventCut && !fEventCut->IsSelected(fESDev, IsCollision())) return;
213
214   if(!fESDfriend){
215     AliError("Failed retrieving ESD friend event");
216     return;
217   }
218   if(HasMCdata() && !fMCev){
219     AliError("Failed retrieving MC event");
220     return;
221   }
222
223   fContainer->Delete();
224   fV0container->Delete();
225   fEventInfo->Delete("");
226   new(fEventInfo)AliTRDeventInfo(fESDev->GetHeader(), const_cast<AliESDRun *>(fESDev->GetESDRun()));
227   
228   Bool_t *trackMap(NULL);
229   AliStack * mStack(NULL);
230   Int_t nTracksMC = HasMCdata() ? fMCev->GetNumberOfTracks() : 0, nTracksESD = fESDev->GetNumberOfTracks();
231   if(HasMCdata()){
232     mStack = fMCev->Stack();
233     if(!mStack){
234       AliError("Failed retrieving MC Stack");
235       return;
236     }
237     trackMap = new Bool_t[nTracksMC];
238     memset(trackMap, 0, sizeof(Bool_t) * nTracksMC);
239   }
240   
241   Double32_t dedx[100]; Int_t nSlices(0);
242   Int_t nTRDout(0), nTRDin(0), nTPC(0), nclsTrklt;
243   AliDebug(2, Form("Entry[%3d] Tracks: ESD[%d] MC[%d]\n", (Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), nTracksESD, nTracksMC));
244   AliESDtrack *esdTrack = NULL;
245   AliESDfriendTrack *esdFriendTrack = NULL;
246   TObject *calObject = NULL;
247   AliTRDtrackV1 *track = NULL;
248   AliTRDseedV1 *tracklet = NULL;
249   AliTRDcluster *cl = NULL;
250   for(Int_t itrk = 0; itrk < nTracksESD; itrk++){
251     esdTrack = fESDev->GetTrack(itrk);
252     AliDebug(3, Form("\n%3d ITS[%d] TPC[%d] TRD[%d]\n", itrk, esdTrack->GetNcls(0), esdTrack->GetNcls(1), esdTrack->GetNcls(2)));
253
254     // Track Selection
255     if(UseLocalTrkSelection()){
256       if(esdTrack->Pt() < fgkPt) continue;
257       if(TMath::Abs(esdTrack->Eta()) > fgkEta) continue;
258       if(esdTrack->GetTPCNcls() < fgkNclTPC) continue;
259       Float_t par[2], cov[3];
260       esdTrack->GetImpactParameters(par, cov);
261       if(IsCollision()){ // cuts on DCA
262         if(TMath::Abs(par[0]) > fgkTrkDCAxy) continue;
263         if(TMath::Abs(par[1]) > fgkTrkDCAz) continue;
264       }
265     }
266     if(fTrackCut && !fTrackCut->IsSelected(esdTrack)) continue;
267
268     if(esdTrack->GetStatus()&AliESDtrack::kTPCout) nTPC++;
269     if(esdTrack->GetStatus()&AliESDtrack::kTRDout) nTRDout++;
270     if(esdTrack->GetStatus()&AliESDtrack::kTRDin) nTRDin++;
271
272     // look at external track param
273     const AliExternalTrackParam *op = esdTrack->GetOuterParam();
274     Double_t xyz[3];
275     if(op){
276       op->GetXYZ(xyz);
277       op->Global2LocalPosition(xyz, op->GetAlpha());
278       AliDebug(3, Form("op @ X[%7.3f]\n", xyz[0]));
279     }
280
281     // read MC info
282     Int_t fPdg = -1;
283     Int_t label = -1; UInt_t alab=UINT_MAX;
284     if(HasMCdata()){
285       label = esdTrack->GetLabel(); 
286       alab = TMath::Abs(label);
287       // register the track
288       if(alab < UInt_t(nTracksMC)){ 
289         trackMap[alab] = kTRUE; 
290       } else { 
291         AliError(Form("MC label[%d] outside scope for Ev[%d] Trk[%d].", label, (Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), itrk));
292         continue; 
293       }
294       AliMCParticle *mcParticle = NULL; 
295       if(!(mcParticle = (AliMCParticle*) fMCev->GetTrack(alab))){
296         AliError(Form("MC particle label[%d] missing for Ev[%d] Trk[%d].", label, (Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), itrk));
297         continue;
298       }
299       fPdg = mcParticle->Particle()->GetPdgCode();
300       Int_t nRefs = mcParticle->GetNumberOfTrackReferences();
301       Int_t iref = 0; AliTrackReference *ref = NULL; 
302       while(iref<nRefs){
303         ref = mcParticle->GetTrackReference(iref);
304         if(ref->LocalX() > fgkTPC) break;
305         iref++;
306       }
307
308       new(fTrackInfo) AliTRDtrackInfo();
309       fTrackInfo->SetPDG(fPdg);
310       fTrackInfo->SetPrimary(mcParticle->Particle()->IsPrimary());
311       Int_t jref = iref;//, kref = 0;
312       while(jref<nRefs){
313         ref = mcParticle->GetTrackReference(jref);
314         if(ref->LocalX() > fgkTOF) break;
315         AliDebug(4, Form("  trackRef[%2d (%2d)] @ %7.3f OK", jref-iref, jref, ref->LocalX()));
316         fTrackInfo->AddTrackRef(ref);
317         jref++;
318       }
319       AliDebug(3, Form("NtrackRefs[%d(%d)]", fTrackInfo->GetNTrackRefs(), nRefs));
320     } else {
321       new (fTrackInfo) AliTRDtrackInfo();
322       fTrackInfo->SetPDG(fPdg);
323     }
324
325     // copy some relevant info to TRD track info
326     fTrackInfo->SetStatus(esdTrack->GetStatus());
327     fTrackInfo->SetTrackId(esdTrack->GetID());
328     Double_t p[AliPID::kSPECIES]; esdTrack->GetTRDpid(p);
329     fTrackInfo->SetESDpid(p);
330     fTrackInfo->SetESDpidQuality(esdTrack->GetTRDntrackletsPID());
331     if(!nSlices) nSlices = esdTrack->GetNumberOfTRDslices();
332     memset(dedx, 0, 100*sizeof(Double32_t));
333     Int_t in(0);
334     for(Int_t il=0; il<AliTRDgeometry::kNlayer; il++)
335       for(Int_t is=0; is<nSlices; is++) 
336         dedx[in++]=esdTrack->GetTRDslice(il, is);
337     for(Int_t il=0; il<AliTRDgeometry::kNlayer; il++) dedx[in++]=esdTrack->GetTRDmomentum(il);
338     fTrackInfo->SetSlices(in, dedx);
339     fTrackInfo->SetLabel(label);
340     fTrackInfo->SetNumberOfClustersRefit(esdTrack->GetNcls(2));
341     // some other Informations which we may wish to store in order to find problematic cases
342     fTrackInfo->SetKinkIndex(esdTrack->GetKinkIndex(0));
343     fTrackInfo->SetTPCncls(static_cast<UShort_t>(esdTrack->GetNcls(1)));
344     nclsTrklt = 0;
345   
346
347     // read REC info
348     esdFriendTrack = fESDfriend->GetTrack(itrk);
349     if(esdFriendTrack){
350       Int_t icalib = 0;
351       while((calObject = esdFriendTrack->GetCalibObject(icalib++))){
352         if(strcmp(calObject->IsA()->GetName(),"AliTRDtrackV1") != 0) continue; // Look for the TRDtrack
353         if(!(track = dynamic_cast<AliTRDtrackV1*>(calObject))) break;
354         AliDebug(4, Form("TRD track OK"));
355         // Set the clusters to unused
356         for(Int_t ipl = 0; ipl < AliTRDgeometry::kNlayer; ipl++){
357           if(!(tracklet = track->GetTracklet(ipl))) continue;
358           tracklet->ResetClusterIter();
359           while((cl = tracklet->NextCluster())) cl->Use(0);
360         }
361         fTrackInfo->SetTrack(track);
362         break;
363       }
364       AliDebug(3, Form("Ntracklets[%d]\n", fTrackInfo->GetNTracklets()));
365     } else AliDebug(3, "No ESD friends");
366     if(op) fTrackInfo->SetOuterParam(op);
367
368     if(DebugLevel() >= 1){
369       AliTRDtrackInfo info(*fTrackInfo);
370       (*DebugStream()) << "trackInfo"
371       << "TrackInfo.=" << &info
372       << "\n";
373       info.Delete("");
374     }
375   
376     fContainer->Add(new AliTRDtrackInfo(*fTrackInfo));
377     fTrackInfo->Delete("");
378   }
379
380   AliDebug(2, Form("%3d Tracks: TPCout[%d] TRDin[%d] TRDout[%d]\n", (Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), nTPC, nTRDin, nTRDout));
381
382 //   AliESDv0 *v0 = NULL;
383 //   for(Int_t iv0=0; iv0<fESD->GetNumberOfV0s(); iv0++){
384 //     if(!(v0 = fESD->GetV0(iv0))) continue;
385 //     fV0container->Add(new AliTRDv0Info(v0));
386 //   }
387
388   // Insert also MC tracks which are passing TRD where the track is not reconstructed
389   if(HasMCdata()){
390     AliDebug(10, "Output of the MC track map:");
391     for(Int_t itk = 0; itk < nTracksMC;  itk++) AliDebug(10, Form("trackMap[%d] = %s", itk, trackMap[itk] == kTRUE ? "TRUE" : "kFALSE"));
392   
393     for(Int_t itk = 0; itk < nTracksMC; itk++){
394       if(trackMap[itk]) continue;
395       AliMCParticle *mcParticle =  (AliMCParticle*) fMCev->GetTrack(TMath::Abs(itk));
396       Int_t fPdg = mcParticle->Particle()->GetPdgCode();
397       Int_t nRefs = mcParticle->GetNumberOfTrackReferences();
398       Int_t iref = 0; AliTrackReference *ref = NULL; 
399       Int_t nRefsTRD = 0;
400       new(fTrackInfo) AliTRDtrackInfo();
401       fTrackInfo->SetPDG(fPdg);
402       while(iref<nRefs){
403         Bool_t kIN(kFALSE);
404         ref = mcParticle->GetTrackReference(iref);
405         if(ref->LocalX() > fgkTPC && ref->LocalX() < fgkTOF){
406           fTrackInfo->AddTrackRef(ref);
407           nRefsTRD++;kIN=kTRUE;
408         }
409         AliDebug(4, Form("  trackRef[%2d] @ x[%7.3f] %s", iref, ref->LocalX(), kIN?"IN":"OUT"));
410         iref++;
411       }
412       if(!nRefsTRD){
413         // In this stage we at least require 1 hit inside TRD. What will be done with this tracks is a task for the 
414         // analysis job
415         fTrackInfo->Delete("");
416         continue;
417       }
418       fTrackInfo->SetPrimary(mcParticle->Particle()->IsPrimary());
419       fTrackInfo->SetLabel(itk);
420       if(DebugLevel() >= 1){
421         AliTRDtrackInfo info(*fTrackInfo);
422         (*DebugStream()) << "trackInfo"
423         << "TrackInfo.=" << &info
424         << "\n";
425         info.Delete("");
426       }
427       AliDebug(3, Form("Add MC track @ label[%d] nTRDrefs[%d].", itk, nRefsTRD));
428       fContainer->Add(new AliTRDtrackInfo(*fTrackInfo));
429       fTrackInfo->Delete("");
430     }
431     delete[] trackMap;
432   }
433   PostData(1, fContainer);
434   PostData(2, fEventInfo);
435   PostData(3, fV0container);
436 }
437