1 /**************************************************************************
\r
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
\r
4 * Author: The ALICE Off-line Project. *
\r
5 * Contributors are mentioned in the code where appropriate. *
\r
7 * Permission to use, copy, modify and distribute this software and its *
\r
8 * documentation strictly for non-commercial purposes is hereby granted *
\r
9 * without fee, provided that the above copyright notice appears in all *
\r
10 * copies and that both the copyright notice and this permission notice *
\r
11 * appear in the supporting documentation. The authors make no claims *
\r
12 * about the suitability of this software for any purpose. It is *
\r
13 * provided "as is" without express or implied warranty. *
\r
14 **************************************************************************/
\r
16 /* $Id: AliTRDinfoGen.cxx 27496 2008-07-22 08:35:45Z cblume $ */
\r
18 ////////////////////////////////////////////////////////////////////////////
\r
20 // Tender wagon for TRD performance/calibration train
\r
22 // In this wagon the information from
\r
24 // - Friends [if available]
\r
25 // - MC [if available]
\r
26 // are grouped into AliTRDtrackInfo objects and fed to worker tasks
\r
29 // Markus Fasel <M.Fasel@gsi.de>
\r
30 // Alexandru Bercuci <A.Bercuci@gsi.de>
\r
32 ////////////////////////////////////////////////////////////////////////////
\r
34 #include <TClonesArray.h>
\r
35 #include <TObjArray.h>
\r
36 #include <TObject.h>
\r
37 #include <TString.h>
\r
45 #include <TParticle.h>
\r
48 #include "AliAnalysisManager.h"
\r
49 #include "AliGeomManager.h"
\r
50 #include "AliCDBManager.h"
\r
51 #include "AliCDBStorage.h"
\r
52 #include "AliCDBEntry.h"
\r
53 #include "AliCDBPath.h"
\r
54 #include "AliESDEvent.h"
\r
55 #include "AliMCEvent.h"
\r
56 #include "AliESDInputHandler.h"
\r
57 #include "AliMCEventHandler.h"
\r
59 #include "AliESDfriend.h"
\r
60 #include "AliESDfriendTrack.h"
\r
61 #include "AliESDHeader.h"
\r
62 #include "AliESDRun.h"
\r
63 #include "AliESDtrack.h"
\r
64 #include "AliESDv0.h"
\r
65 #include "AliESDtrackCuts.h"
\r
66 #include "AliMCParticle.h"
\r
67 #include "AliMultiplicity.h"
\r
69 #include "AliStack.h"
\r
70 #include "AliTrackReference.h"
\r
71 #include "TTreeStream.h"
\r
79 #include "AliTRDReconstructor.h"
\r
80 #include "AliTRDrecoParam.h"
\r
81 #include "AliTRDcalibDB.h"
\r
82 #include "AliTRDtrackerV1.h"
\r
83 #include "AliTRDgeometry.h"
\r
84 #include "AliTRDtrackV1.h"
\r
85 #include "AliTRDseedV1.h"
\r
86 #include "AliTRDcluster.h"
\r
87 #include "AliTRDinfoGen.h"
\r
88 #include "AliTRDpwg1Helper.h"
\r
89 #include "info/AliTRDtrackInfo.h"
\r
90 #include "info/AliTRDeventInfo.h"
\r
91 #include "info/AliTRDv0Info.h"
\r
92 #include "info/AliTRDeventCuts.h"
\r
94 ClassImp(AliTRDinfoGen)
\r
96 const Float_t AliTRDinfoGen::fgkITS = 100.; // to be checked
\r
97 const Float_t AliTRDinfoGen::fgkTPC = 290.;
\r
98 const Float_t AliTRDinfoGen::fgkTRD = 365.;
\r
100 const Float_t AliTRDinfoGen::fgkTrkDCAxy = 3.;
\r
101 const Float_t AliTRDinfoGen::fgkTrkDCAz = 10.;
\r
102 const Int_t AliTRDinfoGen::fgkNclTPC = 70;
\r
103 const Float_t AliTRDinfoGen::fgkPt = 0.2;
\r
104 const Float_t AliTRDinfoGen::fgkEta = 0.9;
\r
105 AliTRDReconstructor* AliTRDinfoGen::fgReconstructor(NULL);
\r
106 AliTRDgeometry* AliTRDinfoGen::fgGeo(NULL);
\r
107 //____________________________________________________________________
\r
108 AliTRDinfoGen::AliTRDinfoGen()
\r
109 :AliAnalysisTaskSE()
\r
115 ,fOCDB("local://$ALICE_ROOT/OCDB")
\r
119 ,fTracksBarrel(NULL)
\r
125 ,fDebugStream(NULL)
\r
128 // Default constructor
\r
130 SetNameTitle("TRDinfoGen", "MC-REC TRD-track list generator");
\r
133 //____________________________________________________________________
\r
134 AliTRDinfoGen::AliTRDinfoGen(char* name)
\r
135 :AliAnalysisTaskSE(name)
\r
141 ,fOCDB("local://$ALICE_ROOT/OCDB")
\r
145 ,fTracksBarrel(NULL)
\r
151 ,fDebugStream(NULL)
\r
154 // Default constructor
\r
156 SetTitle("MC-REC TRD-track list generator");
\r
157 DefineOutput(AliTRDpwg1Helper::kTracksBarrel, TObjArray::Class());
\r
158 DefineOutput(AliTRDpwg1Helper::kTracksSA, TObjArray::Class());
\r
159 DefineOutput(AliTRDpwg1Helper::kTracksKink, TObjArray::Class());
\r
160 DefineOutput(AliTRDpwg1Helper::kEventInfo, AliTRDeventInfo::Class());
\r
161 DefineOutput(AliTRDpwg1Helper::kV0List, TObjArray::Class());
\r
162 DefineOutput(AliTRDpwg1Helper::kMonitor, TObjArray::Class()); // histogram list
\r
165 //____________________________________________________________________
\r
166 AliTRDinfoGen::~AliTRDinfoGen()
\r
169 if(fgGeo) delete fgGeo;
\r
170 if(fgReconstructor) delete fgReconstructor;
\r
171 if(fDebugStream) delete fDebugStream;
\r
172 if(fV0Cut) delete fV0Cut;
\r
173 if(fTrackCut) delete fTrackCut;
\r
174 if(fEventCut) delete fEventCut;
\r
175 if(fTrackInfo) delete fTrackInfo; fTrackInfo = NULL;
\r
176 if(fEventInfo) delete fEventInfo; fEventInfo = NULL;
\r
177 if(fV0Info) delete fV0Info; fV0Info = NULL;
\r
178 if(fTracksBarrel){
\r
179 fTracksBarrel->Delete(); delete fTracksBarrel;
\r
180 fTracksBarrel = NULL;
\r
183 fTracksSA->Delete(); delete fTracksSA;
\r
187 fTracksKink->Delete(); delete fTracksKink;
\r
188 fTracksKink = NULL;
\r
191 fV0List->Delete();
\r
195 if(fContainer && !(AliAnalysisManager::GetAnalysisManager() && AliAnalysisManager::GetAnalysisManager()->IsProofMode())){
\r
196 fContainer->Delete();
\r
202 //____________________________________________________________________
\r
203 Bool_t AliTRDinfoGen::GetRefFigure(Int_t)
\r
205 // General graphs for PWG1/TRD train
\r
207 AliWarning("Please provide a canvas to draw results.");
\r
210 fContainer->At(kStatTrk)->Draw("bar");
\r
214 //____________________________________________________________________
\r
215 void AliTRDinfoGen::UserCreateOutputObjects()
\r
218 // Create Output Containers (TObjectArray containing 1D histograms)
\r
221 fTrackInfo = new AliTRDtrackInfo();
\r
222 fEventInfo = new AliTRDeventInfo();
\r
223 fV0Info = new AliTRDv0Info();
\r
224 fTracksBarrel = new TObjArray(200); fTracksBarrel->SetOwner(kTRUE);
\r
225 fTracksSA = new TObjArray(20); fTracksSA->SetOwner(kTRUE);
\r
226 fTracksKink = new TObjArray(20); fTracksKink->SetOwner(kTRUE);
\r
227 fV0List = new TObjArray(10); fV0List->SetOwner(kTRUE);
\r
229 // define general monitor
\r
230 fContainer = new TObjArray(kNclasses); fContainer->SetOwner(kTRUE);
\r
231 TH1 *h=new TH1I("hStat", "Run statistics;Observable;Entries", Int_t(kNObjects), -0.5, Float_t(kNObjects)-0.5);
\r
232 TAxis *ax(h->GetXaxis());
\r
233 ax->SetBinLabel(Int_t(kTracksESD) + 1, "ESD");
\r
234 ax->SetBinLabel(Int_t(kTracksMC) + 1, "MC");
\r
235 ax->SetBinLabel(Int_t(kV0) + 1, "V0");
\r
236 ax->SetBinLabel(Int_t(kTPC) + 1, "TPC");
\r
237 ax->SetBinLabel(Int_t(kTRDin) + 1, "TRDin");
\r
238 ax->SetBinLabel(Int_t(kTRDout) + 1, "TRDout");
\r
239 ax->SetBinLabel(Int_t(kBarrel) + 1, "Barrel");
\r
240 ax->SetBinLabel(Int_t(kBarrelMC) + 1, "BarrelMC");
\r
241 ax->SetBinLabel(Int_t(kSA) + 1, "SA");
\r
242 ax->SetBinLabel(Int_t(kSAMC) + 1, "SAMC");
\r
243 ax->SetBinLabel(Int_t(kKink) + 1, "Kink");
\r
244 ax->SetBinLabel(Int_t(kKinkMC) + 1, "KinkMC");
\r
245 ax->SetBinLabel(Int_t(kBarrelFriend) + 1, "BFriend");
\r
246 ax->SetBinLabel(Int_t(kSAFriend) + 1, "SFriend");
\r
247 fContainer->AddAt(h, kStatTrk);
\r
248 h=new TH1I("hEv", "Run statistics;Event Class;Entries", 4, -0.5, 3.5);
\r
249 ax = h->GetXaxis();
\r
250 ax->SetBinLabel(1, "Low");
\r
251 ax->SetBinLabel(2, "High");
\r
252 ax->SetBinLabel(3, "Cosmic");
\r
253 ax->SetBinLabel(4, "Calib");
\r
254 fContainer->AddAt(h, kEvType);
\r
255 TH2I* h2=new TH2I("hBC", "Bunch Cross statistics;Fill Bunch;TOF BC;Entries", 3500, -0.5, 3499.5, 31, -10.5, 20.5);
\r
256 fContainer->AddAt(h2, kBC);
\r
257 h=new TH1I("hTriggers", "Triggers statistics;;Entries", 21, -0.5, 20.5);
\r
258 fContainer->AddAt(h, kTrigger);
\r
259 PostData(AliTRDpwg1Helper::kMonitor, fContainer);
\r
262 //____________________________________________________________________
\r
263 Bool_t AliTRDinfoGen::Load(const Char_t *file, const Char_t *dir, const Char_t *name)
\r
265 // Load data from performance file
\r
267 if(!TFile::Open(file)){
\r
268 AliWarning(Form("Couldn't open file %s.", file));
\r
272 if(!gFile->cd(dir)){
\r
273 AliWarning(Form("Couldn't cd to %s in %s.", dir, file));
\r
277 TObjArray *o(NULL);
\r
278 const Char_t *tn=(name ? name : GetName());
\r
279 if(!(o = (TObjArray*)gDirectory->Get(tn))){
\r
280 AliWarning(Form("Missing histogram container %s.", tn));
\r
283 fContainer = (TObjArray*)o->Clone(GetName());
\r
288 //____________________________________________________________________
\r
289 void AliTRDinfoGen::UserExec(Option_t *){
\r
291 // Run the Analysis
\r
294 fTracksBarrel->Delete();
\r
295 fTracksSA->Delete();
\r
296 fTracksKink->Delete();
\r
298 fEventInfo->Delete("");
\r
300 fESDev = dynamic_cast<AliESDEvent*>(InputEvent());
\r
302 AliError("Failed retrieving ESD event");
\r
306 // This part may conflict with other detectors !!
\r
308 AliCDBEntry* obj(NULL);
\r
309 AliCDBManager* ocdb = AliCDBManager::Instance();
\r
310 if(ocdb->IsDefaultStorageSet()){
\r
311 AliInfo("OCDB :: initialized.");
\r
313 AliInfo("OCDB :: initializing locally ...");
\r
314 // prepare OCDB access
\r
315 ocdb->SetDefaultStorage(fOCDB.Data());
\r
316 ocdb->SetRun(fESDev->GetRunNumber());
\r
317 // create geo manager
\r
318 if(!(obj = ocdb->Get(AliCDBPath("GRP", "Geometry", "Data")))){
\r
319 AliError("GEOMETRY failed initialization.");
\r
321 AliGeomManager::SetGeometry((TGeoManager*)obj->GetObject());
\r
322 AliGeomManager::GetNalignable("TRD");
\r
323 AliGeomManager::ApplyAlignObjsFromCDB("TRD");
\r
325 //init magnetic field
\r
326 if(!TGeoGlobalMagField::Instance()->IsLocked() &&
\r
327 !fESDev->InitMagneticField()){
\r
328 AliError("MAGNETIC FIELD failed initialization.");
\r
330 // set no of time bins
\r
331 AliTRDtrackerV1::SetNTimeBins(AliTRDcalibDB::Instance()->GetNumberOfTimeBinsDCS());
\r
333 AliInfo(Form("OCDB : Loc[%s] Run[%d] TB[%d]", ocdb->GetDefaultStorage()->GetURI().Data(), ocdb->GetRun(), AliTRDtrackerV1::GetNTimeBins()));
\r
335 // load misalignment
\r
336 fgGeo = new AliTRDgeometry;
\r
337 fgGeo->CreateClusterMatrixArray();
\r
338 // load reco param list from OCDB
\r
339 AliInfo("Initializing TRD reco params ...");
\r
340 fgReconstructor = new AliTRDReconstructor();
\r
341 if(!(obj = ocdb->Get(AliCDBPath("TRD", "Calib", "RecoParam")))){
\r
342 AliError("RECO PARAM failed initialization.");
\r
344 obj->PrintMetaData();
\r
345 fRecos = (TObjArray*)obj->GetObject();
\r
350 AliDebug(2, "+++++++++++++++++++++++++++++++++++++++");
\r
351 AliDebug(2, Form("Analyse Ev[%d]", fESDev->GetEventNumberInFile()));
\r
353 // set reco param valid for this event
\r
354 TH1 *h = (TH1I*)fContainer->At(kEvType);
\r
356 fgReconstructor->SetRecoParam(AliTRDrecoParam::GetLowFluxParam());
\r
359 for(Int_t ireco(0); ireco<fRecos->GetEntriesFast(); ireco++){
\r
360 AliTRDrecoParam *reco((AliTRDrecoParam*)fRecos->At(ireco));
\r
361 Int_t es(reco->GetEventSpecie());
\r
362 if(!(es&fESDev->GetEventSpecie())) continue;
\r
363 fgReconstructor->SetRecoParam(reco);
\r
364 if(AliLog::GetDebugLevel("PWG1/TRD", "AliTRDinfoGen")>5) reco->Dump();
\r
366 if(es&AliRecoParam::kLowMult){ s="LowMult"; h->Fill(0);}
\r
367 else if(es&AliRecoParam::kHighMult){ s="HighMult"; h->Fill(1);}
\r
368 else if(es&AliRecoParam::kCosmic){ s="Cosmic"; h->Fill(2);}
\r
369 else if(es&AliRecoParam::kCalib){ s="Calib"; h->Fill(3);}
\r
371 AliDebug(2, Form(" Using reco param \"%s\" for event %d.", s.Data(), fESDev->GetEventNumberInFile()));
\r
376 // link MC if available
\r
380 h = (TH1I*)fContainer->At(kTrigger);
\r
381 TAxis *ax(h->GetXaxis());
\r
382 TObjArray *evTriggers = fESDev->GetFiredTriggerClasses().Tokenize(" ");
\r
383 for(Int_t iet(evTriggers->GetEntriesFast()); iet--;){
\r
385 for(; ix<=ax->GetNbins(); ix++){
\r
386 if(!Int_t(h->GetBinContent(ix))){
\r
387 ax->SetBinLabel(ix, (*evTriggers)[iet]->GetName());
\r
390 if(strcmp((*evTriggers)[iet]->GetName(), ax->GetBinLabel(ix))==0) break;
\r
392 h->AddBinContent(ix);
\r
395 // event selection based on vertex cuts and trigger
\r
396 if(UseLocalEvSelection() && !fEventCut->IsSelected(fESDev, IsCollision())) return;
\r
399 AliError("Failed retrieving ESD friend event");
\r
402 if(HasMCdata() && !fMCev){
\r
403 AliError("Failed retrieving MC event");
\r
407 new(fEventInfo)AliTRDeventInfo(fESDev->GetHeader(), const_cast<AliESDRun *>(fESDev->GetESDRun()));
\r
408 // Determine centrality
\r
409 // Author: Ionut Arsene <I.C.Arsene@gsi.de>
\r
410 Int_t centralityBin = -1;
\r
411 AliDebug(2, Form(" Beam Type: %s", fESDev->GetESDRun()->GetBeamType()));
\r
412 TString beamtype = fESDev->GetESDRun()->GetBeamType();
\r
413 if(beamtype.Contains("Pb-Pb") || beamtype.Contains("A-A")){
\r
415 const AliMultiplicity *mult = fESDev->GetMultiplicity();
\r
416 Double_t zdcNeutronEnergy = fESDev->GetZDCN1Energy()+fESDev->GetZDCN2Energy();
\r
417 Double_t itsNTracklets = mult->GetNumberOfTracklets();
\r
418 Double_t centralitySlopes[6] = {0.0, 4.0, 8.0, 20.0, 50.0, 1000000.};
\r
419 AliDebug(1, Form("zdcNeutronEnergy: %f, itsNTracklets: %f\n", zdcNeutronEnergy, itsNTracklets));
\r
420 for(Int_t iCent=1; iCent<=5; ++iCent) {
\r
421 if(zdcNeutronEnergy>centralitySlopes[iCent-1]*itsNTracklets && zdcNeutronEnergy<centralitySlopes[iCent]*itsNTracklets)
\r
422 centralityBin=iCent - 1;
\r
424 AliDebug(2, Form(" Centrality Class: %d", centralityBin));
\r
426 fEventInfo->SetCentrality(centralityBin);
\r
427 UShort_t evBC(fESDev->GetBunchCrossNumber());
\r
429 Bool_t *trackMap(NULL);
\r
430 AliStack * mStack(NULL);
\r
431 Int_t nTracksMC = HasMCdata() ? fMCev->GetNumberOfTracks() : 0, nTracksESD = fESDev->GetNumberOfTracks();
\r
433 mStack = fMCev->Stack();
\r
435 AliError("Failed retrieving MC Stack");
\r
438 trackMap = new Bool_t[nTracksMC];
\r
439 memset(trackMap, 0, sizeof(Bool_t) * nTracksMC);
\r
442 Double32_t dedx[100]; Int_t nSlices(0);
\r
443 Int_t nTRDout(0), nTRDin(0), nTPC(0)
\r
445 ,nBarrel(0), nSA(0), nKink(0)
\r
446 ,nBarrelFriend(0), nSAFriend(0)
\r
447 ,nBarrelMC(0), nSAMC(0), nKinkMC(0);
\r
448 AliESDtrack *esdTrack = NULL;
\r
449 AliESDfriendTrack *esdFriendTrack = NULL;
\r
450 TObject *calObject = NULL;
\r
451 AliTRDtrackV1 *track = NULL;
\r
452 AliTRDseedV1 *tracklet = NULL;
\r
453 AliTRDcluster *cl = NULL;
\r
456 // LOOP 0 - over ESD v0s
\r
457 Float_t bField(fESDev->GetMagneticField());
\r
458 AliESDv0 *v0(NULL);
\r
459 Int_t v0pid[AliPID::kSPECIES];
\r
460 for(Int_t iv0(0); iv0<fESDev->GetNumberOfV0s(); iv0++){
\r
461 if(!(v0 = fESDev->GetV0(iv0))) continue;
\r
463 if(fV0Cut) new(fV0Info) AliTRDv0Info(*fV0Cut);
\r
464 else new(fV0Info) AliTRDv0Info();
\r
465 fV0Info->SetMagField(bField);
\r
466 fV0Info->SetV0tracks(fESDev->GetTrack(v0->GetPindex()), fESDev->GetTrack(v0->GetNindex()));
\r
467 fV0Info->SetV0Info(v0);
\r
468 fV0List->Add(new AliTRDv0Info(*fV0Info));// kFOUND=kFALSE;
\r
472 // LOOP 1 - over ESD tracks
\r
473 AliTRDv0Info *v0info=NULL;
\r
474 for(Int_t itrk = 0; itrk < nTracksESD; itrk++){
\r
475 new(fTrackInfo) AliTRDtrackInfo();
\r
476 esdTrack = fESDev->GetTrack(itrk);
\r
477 AliDebug(3, Form("\n%3d ITS[%d] TPC[%d] TRD[%d] TOF-BC[%d]\n", itrk, esdTrack->GetNcls(0), esdTrack->GetNcls(1), esdTrack->GetNcls(2), esdTrack->GetTOFBunchCrossing()));
\r
479 if(esdTrack->GetStatus()&AliESDtrack::kTPCout) nTPC++;
\r
480 if(esdTrack->GetStatus()&AliESDtrack::kTRDout) nTRDout++;
\r
481 if(esdTrack->GetStatus()&AliESDtrack::kTRDin) nTRDin++;
\r
483 // look at external track param
\r
484 const AliExternalTrackParam *op = esdTrack->GetOuterParam();
\r
488 op->Global2LocalPosition(xyz, op->GetAlpha());
\r
489 AliDebug(3, Form("op @ X[%7.3f]\n", xyz[0]));
\r
494 Int_t label = -1; UInt_t alab=UINT_MAX;
\r
496 label = esdTrack->GetLabel();
\r
497 alab = TMath::Abs(label);
\r
498 // register the track
\r
499 if(alab < UInt_t(nTracksMC)){
\r
500 trackMap[alab] = kTRUE;
\r
502 AliError(Form("MC label[%d] outside scope for Ev[%d] Trk[%d].", label, (Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), itrk));
\r
505 AliMCParticle *mcParticle = NULL;
\r
506 if(!(mcParticle = (AliMCParticle*) fMCev->GetTrack(alab))){
\r
507 AliError(Form("MC particle label[%d] missing for Ev[%d] Trk[%d].", label, (Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), itrk));
\r
510 fPdg = mcParticle->Particle()->GetPdgCode();
\r
511 Int_t nRefs = mcParticle->GetNumberOfTrackReferences();
\r
512 Int_t iref = 0; AliTrackReference *ref = NULL;
\r
514 ref = mcParticle->GetTrackReference(iref);
\r
515 if(ref->LocalX() > fgkTPC) break;
\r
519 fTrackInfo->SetMC();
\r
520 fTrackInfo->SetPDG(fPdg);
\r
521 fTrackInfo->SetPrimary(mcParticle->Particle()->IsPrimary());
\r
522 fTrackInfo->SetLabel(label);
\r
523 Int_t jref = iref;//, kref = 0;
\r
525 ref = mcParticle->GetTrackReference(jref);
\r
526 if(ref->LocalX() > fgkTRD) break;
\r
527 AliDebug(4, Form(" trackRef[%2d (%2d)] @ %7.3f OK", jref-iref, jref, ref->LocalX()));
\r
528 fTrackInfo->AddTrackRef(ref);
\r
531 AliDebug(3, Form("NtrackRefs[%d(%d)]", fTrackInfo->GetNTrackRefs(), nRefs));
\r
534 // copy some relevant info to TRD track info
\r
535 fTrackInfo->SetStatus(esdTrack->GetStatus());
\r
536 fTrackInfo->SetTrackId(esdTrack->GetID());
\r
537 Double_t p[AliPID::kSPECIES]; esdTrack->GetTRDpid(p);
\r
538 fTrackInfo->SetESDpid(p);
\r
539 fTrackInfo->SetESDpidQuality(esdTrack->GetTRDntrackletsPID());
\r
540 if(!nSlices) nSlices = esdTrack->GetNumberOfTRDslices();
\r
541 memset(dedx, 0, 100*sizeof(Double32_t));
\r
543 for(Int_t il=0; il<AliTRDgeometry::kNlayer; il++)
\r
544 for(Int_t is=0; is<nSlices; is++)
\r
545 dedx[in++]=esdTrack->GetTRDslice(il, is);
\r
546 for(Int_t il=0; il<AliTRDgeometry::kNlayer; il++) dedx[in++]=esdTrack->GetTRDmomentum(il);
\r
547 fTrackInfo->SetSlices(in, dedx);
\r
548 fTrackInfo->SetNumberOfClustersRefit(esdTrack->GetNcls(2));
\r
549 // some other Informations which we may wish to store in order to find problematic cases
\r
550 fTrackInfo->SetKinkIndex(esdTrack->GetKinkIndex(0));
\r
551 fTrackInfo->SetTPCncls(static_cast<UShort_t>(esdTrack->GetNcls(1)));
\r
552 fTrackInfo->SetTOFbc(esdTrack->GetTOFBunchCrossing()==AliVTrack::kTOFBCNA?0:esdTrack->GetTOFBunchCrossing());
\r
556 for(Int_t iv(0); iv<fV0List->GetEntriesFast(); iv++){
\r
557 if(!(v0info = (AliTRDv0Info*)fV0List->At(iv))) continue;
\r
558 if(!v0info->GetV0Daughter(1) && !v0info->GetV0Daughter(-1)) continue;
\r
559 if(!v0info->HasTrack(fTrackInfo)) continue;
\r
560 memset(v0pid, 0, AliPID::kSPECIES*sizeof(Int_t));
\r
561 fTrackInfo->SetV0();
\r
562 for(Int_t is=AliPID::kSPECIES; is--;){v0pid[is] = v0info->GetPID(is, fTrackInfo);}
\r
563 fTrackInfo->SetV0pid(v0pid);
\r
564 fTrackInfo->SetV0();
\r
565 //const AliTRDtrackInfo::AliESDinfo *ei = fTrackInfo->GetESDinfo();
\r
570 esdFriendTrack = (fESDfriend->GetNumberOfTracks() > itrk) ? fESDfriend->GetTrack(itrk): NULL;
\r
572 if(esdFriendTrack){
\r
573 fTrackInfo->SetTPCoutParam(esdFriendTrack->GetTPCOut());
\r
575 while((calObject = esdFriendTrack->GetCalibObject(icalib++))){
\r
576 if(strcmp(calObject->IsA()->GetName(),"AliTRDtrackV1") != 0) continue; // Look for the TRDtrack
\r
577 if(!(track = dynamic_cast<AliTRDtrackV1*>(calObject))) break;
\r
578 AliDebug(4, Form("TRD track OK"));
\r
579 // Set the clusters to unused
\r
580 for(Int_t ipl = 0; ipl < AliTRDgeometry::kNlayer; ipl++){
\r
581 if(!(tracklet = track->GetTracklet(ipl))) continue;
\r
582 tracklet->ResetClusterIter();
\r
583 while((cl = tracklet->NextCluster())) cl->Use(0);
\r
585 fTrackInfo->SetTrack(track);
\r
586 ((TH2I*)fContainer->At(kBC))->Fill(evBC, esdTrack->GetTOFBunchCrossing());
\r
589 AliDebug(3, Form("Ntracklets[%d]\n", fTrackInfo->GetNTracklets()));
\r
590 } else AliDebug(3, "No ESD friends");
\r
591 if(op) fTrackInfo->SetOuterParam(op);
\r
593 if(DebugLevel() >= 1){
\r
594 AliTRDtrackInfo info(*fTrackInfo);
\r
595 (*DebugStream()) << "trackInfo"
\r
596 << "TrackInfo.=" << &info
\r
601 ULong_t status(esdTrack->GetStatus());
\r
602 if((status&AliESDtrack::kTPCout)){
\r
603 if(!esdTrack->GetKinkIndex(0)){ // Barrel Track Selection
\r
604 Bool_t selected(kTRUE);
\r
605 if(UseLocalTrkSelection()){
\r
606 if(esdTrack->Pt() < fgkPt){
\r
607 AliDebug(3, Form("Reject Trk[%3d] Ev[%4d] Pt[%5.2f]", itrk, fESDev->GetEventNumberInFile(), esdTrack->Pt()));
\r
610 if(selected && TMath::Abs(esdTrack->Eta()) > fgkEta){
\r
611 AliDebug(3, Form("Reject Trk[%3d] Ev[%4d] Eta[%5.2f]", itrk, fESDev->GetEventNumberInFile(), TMath::Abs(esdTrack->Eta())));
\r
614 if(selected && esdTrack->GetTPCNcls() < fgkNclTPC){
\r
615 AliDebug(3, Form("Reject Trk[%3d] Ev[%4d] NclTPC[%d]", itrk, fESDev->GetEventNumberInFile(), esdTrack->GetTPCNcls()));
\r
618 Float_t par[2], cov[3];
\r
619 esdTrack->GetImpactParameters(par, cov);
\r
620 if(IsCollision()){ // cuts on DCA
\r
621 if(selected && TMath::Abs(par[0]) > fgkTrkDCAxy){
\r
622 AliDebug(3, Form("Reject Trk[%3d] Ev[%4d] DCAxy[%f]", itrk, fESDev->GetEventNumberInFile(), TMath::Abs(par[0])));
\r
625 if(selected && TMath::Abs(par[1]) > fgkTrkDCAz){
\r
626 AliDebug(3, Form("Reject Trk[%3d] Ev[%4d] DCAz[%f]", itrk, fESDev->GetEventNumberInFile(), TMath::Abs(par[1])));
\r
629 } else if(selected && fMCev && !fMCev->IsPhysicalPrimary(alab)){;
\r
630 AliDebug(3, Form("Reject Trk[%3d] Ev[%4d] Primary", itrk, fESDev->GetEventNumberInFile()));
\r
634 if(fTrackCut && !fTrackCut->IsSelected(esdTrack)) selected = kFALSE;
\r
636 fTracksBarrel->Add(new AliTRDtrackInfo(*fTrackInfo));
\r
638 if(fTrackInfo->GetTrack())
\r
642 fTracksKink->Add(new AliTRDtrackInfo(*fTrackInfo));
\r
645 } else if((status&AliESDtrack::kTRDout) && !(status&AliESDtrack::kTRDin)){
\r
646 fTracksSA->Add(new AliTRDtrackInfo(*fTrackInfo));
\r
648 if(fTrackInfo->GetTrack())
\r
651 fTrackInfo->Delete("");
\r
654 // LOOP 2 - over MC tracks which are passing TRD where the track is not reconstructed
\r
656 AliDebug(10, "Output of the MC track map:");
\r
657 for(Int_t itk = 0; itk < nTracksMC; itk++) AliDebug(10, Form("trackMap[%d] = %s", itk, trackMap[itk] == kTRUE ? "TRUE" : "kFALSE"));
\r
659 for(Int_t itk = 0; itk < nTracksMC; itk++){
\r
660 if(trackMap[itk]) continue;
\r
661 AliMCParticle *mcParticle = (AliMCParticle*) fMCev->GetTrack(TMath::Abs(itk));
\r
662 Int_t fPdg = mcParticle->Particle()->GetPdgCode();
\r
663 Int_t nRefs = mcParticle->GetNumberOfTrackReferences();
\r
664 Int_t iref = 0; AliTrackReference *ref = NULL;
\r
665 Int_t nRefsTRD = 0;
\r
666 new(fTrackInfo) AliTRDtrackInfo();
\r
667 fTrackInfo->SetMC();
\r
668 fTrackInfo->SetPDG(fPdg);
\r
669 while(iref<nRefs){ // count TRD TR
\r
670 Bool_t kIN(kFALSE);
\r
671 ref = mcParticle->GetTrackReference(iref);
\r
672 if(ref->LocalX() > fgkTPC && ref->LocalX() < fgkTRD){
\r
673 fTrackInfo->AddTrackRef(ref);
\r
674 nRefsTRD++;kIN=kTRUE;
\r
676 AliDebug(4, Form(" trackRef[%2d] @ x[%7.3f] %s", iref, ref->LocalX(), kIN?"IN":"OUT"));
\r
680 // In this stage we at least require 1 hit inside TRD. What will be done with this tracks is a task for the
\r
682 fTrackInfo->Delete("");
\r
685 fTrackInfo->SetPrimary(mcParticle->Particle()->IsPrimary());
\r
686 fTrackInfo->SetLabel(itk);
\r
687 if(DebugLevel() >= 1){
\r
688 AliTRDtrackInfo info(*fTrackInfo);
\r
689 (*DebugStream()) << "trackInfo"
\r
690 << "TrackInfo.=" << &info
\r
694 AliDebug(3, Form("Add MC track @ label[%d] nTRDrefs[%d].", itk, nRefsTRD));
\r
695 // check where the track starts
\r
696 ref = mcParticle->GetTrackReference(0);
\r
697 if(ref->LocalX() < fgkITS){
\r
698 fTracksBarrel->Add(new AliTRDtrackInfo(*fTrackInfo));
\r
700 } else if(ref->LocalX() < fgkTPC) {
\r
701 fTracksKink->Add(new AliTRDtrackInfo(*fTrackInfo));
\r
703 } else if(nRefsTRD>6){
\r
704 fTracksSA->Add(new AliTRDtrackInfo(*fTrackInfo));
\r
707 fTrackInfo->Delete("");
\r
712 "\nEv[%3d] Tracks: ESD[%d] MC[%d] V0[%d]\n"
\r
713 " TPCout[%d] TRDin[%d] TRDout[%d]\n"
\r
714 " Barrel[%3d+%3d=%3d] SA[%2d+%2d=%2d] Kink[%2d+%2d=%2d]"
\r
715 ,(Int_t)AliAnalysisManager::GetAnalysisManager()->GetCurrentEntry(), nTracksESD, nTracksMC, fV0List->GetEntries()
\r
716 , nTPC, nTRDin, nTRDout
\r
717 ,nBarrel, nBarrelMC, fTracksBarrel->GetEntries()
\r
718 ,nSA, nSAMC, fTracksSA->GetEntries()
\r
719 ,nKink, nKinkMC, fTracksKink->GetEntries()
\r
721 // save track statistics
\r
722 h = (TH1I*)fContainer->At(kStatTrk);
\r
723 h->Fill(Float_t(kTracksESD), nTracksESD);
\r
724 h->Fill(Float_t(kTracksMC), nTracksMC);
\r
725 h->Fill(Float_t(kV0), fV0List->GetEntries());
\r
726 h->Fill(Float_t(kTPC), nTPC);
\r
727 h->Fill(Float_t(kTRDin), nTRDin);
\r
728 h->Fill(Float_t(kTRDout), nTRDout);
\r
729 h->Fill(Float_t(kBarrel), nBarrel);
\r
730 h->Fill(Float_t(kBarrelMC), nBarrelMC);
\r
731 h->Fill(Float_t(kSA), nSA);
\r
732 h->Fill(Float_t(kSAMC), nSAMC);
\r
733 h->Fill(Float_t(kKink), nKink);
\r
734 h->Fill(Float_t(kKinkMC), nKinkMC);
\r
735 h->Fill(Float_t(kBarrelFriend), nBarrelFriend);
\r
736 h->Fill(Float_t(kSAFriend), nSAFriend);
\r
738 PostData(AliTRDpwg1Helper::kTracksBarrel, fTracksBarrel);
\r
739 PostData(AliTRDpwg1Helper::kTracksSA, fTracksSA);
\r
740 PostData(AliTRDpwg1Helper::kTracksKink, fTracksKink);
\r
741 PostData(AliTRDpwg1Helper::kEventInfo, fEventInfo);
\r
742 PostData(AliTRDpwg1Helper::kV0List, fV0List);
\r
745 //____________________________________________________________________
\r
746 void AliTRDinfoGen::SetLocalEvSelection(const AliTRDeventCuts &ec)
\r
748 // Set event cuts from outside
\r
749 if(!fEventCut) fEventCut = new AliTRDeventCuts(ec);
\r
750 else new(fEventCut) AliTRDeventCuts(ec);
\r
751 fEventCut->Print();
\r
754 //____________________________________________________________________
\r
755 void AliTRDinfoGen::SetLocalV0Selection(const AliTRDv0Info &v0)
\r
757 // Set V0 cuts from outside
\r
759 if(!fV0Cut) fV0Cut = new AliTRDv0Info(v0);
\r
760 else new(fV0Cut) AliTRDv0Info(v0);
\r
764 //____________________________________________________________________
\r
765 TTreeSRedirector* AliTRDinfoGen::DebugStream()
\r
767 // Manage debug stream for task
\r
769 TDirectory *savedir = gDirectory;
\r
770 fDebugStream = new TTreeSRedirector("TRD.DebugInfoGen.root");
\r
773 return fDebugStream;
\r
776 //____________________________________________________________________
\r
777 void AliTRDinfoGen::Terminate(Option_t* /*option*/)
\r
779 // Process run information
\r
781 if(!(fContainer = dynamic_cast<TObjArray *>(GetOutputData(AliTRDpwg1Helper::kMonitor)))) return;
\r
782 AliInfo(Form("fContainer(%p)", (void*)fContainer));
\r
783 if(UseLocalEvSelection()){
\r
784 TH1 *h1 = (TH1*)fContainer->At(kTrigger); TAxis *ax(h1->GetXaxis());
\r
785 AliInfo(Form("h1(%p)", (void*)h1));
\r
786 for(Int_t ix(1); ix<=ax->GetNbins(); ix++){
\r
787 if(fEventCut->CheckTrigger(ax->GetBinLabel(ix))) ax->SetBinLabel(ix, Form("#color[2]{%s}", ax->GetBinLabel(ix)));
\r