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