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