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