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