]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG3/muondep/AliAnalysisTaskMuonFakes.cxx
New task to study fake tracks in simulation: need to access the recoParam used during...
[u/mrichter/AliRoot.git] / PWG3 / muondep / AliAnalysisTaskMuonFakes.cxx
CommitLineData
fa97374c 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// ROOT includes
17#include <TObjArray.h>
18#include <TClonesArray.h>
19#include <TH1F.h>
20#include <TH2F.h>
21#include <TCanvas.h>
22
23// STEER includes
24#include "AliLog.h"
25#include "AliESDEvent.h"
26#include "AliESDMuonTrack.h"
27#include "AliESDInputHandler.h"
28#include "AliMCEventHandler.h"
29
30// ANALYSIS includes
31#include "AliAnalysisManager.h"
32#include "AliCDBManager.h"
33
34// MUON includes
35#include "AliMUONCDB.h"
36#include "AliMUONRecoParam.h"
37#include "AliMUONRecoCheck.h"
38#include "AliMUONVCluster.h"
39#include "AliMUONVTrackStore.h"
40#include "AliMUONTrack.h"
41#include "AliMUONTrackParam.h"
42#include "AliMUONESDInterface.h"
43
44#include "AliAnalysisTaskMuonFakes.h"
45#include "AliCounterCollection.h"
46
47ClassImp(AliAnalysisTaskMuonFakes)
48
49//________________________________________________________________________
50AliAnalysisTaskMuonFakes::AliAnalysisTaskMuonFakes() :
51AliAnalysisTaskSE(),
52fList(0x0),
53fCanvases(0x0),
54fTrackCounters(0x0),
55fFakeTrackCounters(0x0),
56fMatchedTrackCounters(0x0),
57fEventCounters(0x0),
58fCurrentFileName(""),
59fRequestedStationMask(0),
60fRequest2ChInSameSt45(kFALSE),
61fSigmaCut(-1.),
62fUseLabel(kFALSE),
63fRecoParamLocation("alien://folder=/alice/simulation/2008/v4-15-Release/Full")
64{
65 /// Default constructor.
66}
67
68//________________________________________________________________________
69AliAnalysisTaskMuonFakes::AliAnalysisTaskMuonFakes(const char *name) :
70AliAnalysisTaskSE(name),
71fList(0x0),
72fCanvases(0x0),
73fTrackCounters(0x0),
74fFakeTrackCounters(0x0),
75fMatchedTrackCounters(0x0),
76fEventCounters(0x0),
77fCurrentFileName(""),
78fRequestedStationMask(0),
79fRequest2ChInSameSt45(kFALSE),
80fSigmaCut(-1.),
81fUseLabel(kFALSE),
82fRecoParamLocation("alien://folder=/alice/simulation/2008/v4-15-Release/Full")
83{
84 /// Constructor.
85 // Output slot #1 writes into a TObjArray container
86 DefineOutput(1,TObjArray::Class());
87 // Output slot #2 writes into an AliCounterCollection container
88 DefineOutput(2,AliCounterCollection::Class());
89 // Output slot #3 writes into an AliCounterCollection container
90 DefineOutput(3,AliCounterCollection::Class());
91 // Output slot #4 writes into an AliCounterCollection container
92 DefineOutput(4,AliCounterCollection::Class());
93 // Output slot #5 writes into an AliCounterCollection container
94 DefineOutput(5,AliCounterCollection::Class());
95}
96
97//________________________________________________________________________
98AliAnalysisTaskMuonFakes::~AliAnalysisTaskMuonFakes()
99{
100 /// Destructor.
101 if (!AliAnalysisManager::GetAnalysisManager()->IsProofMode()) {
102 delete fList;
103 delete fTrackCounters;
104 delete fFakeTrackCounters;
105 delete fMatchedTrackCounters;
106 delete fEventCounters;
107 }
108 delete fCanvases;
109}
110
111//___________________________________________________________________________
112void AliAnalysisTaskMuonFakes::UserCreateOutputObjects()
113{
114 /// Create histograms and counters.
115
116 fList = new TObjArray(100);
117 fList->SetOwner();
118
119 // number of tracks
120 TH1F *hNumberOfTracks = new TH1F("hNumberOfTracks", "nb of tracks /evt", 21, -0.5, 20.5);
121 fList->AddAtAndExpand(hNumberOfTracks, kNumberOfTracks);
122 TH1F *hNumberOfAdditionalTracks = new TH1F("hNumberOfAdditionalTracks", "nb of fake - nb of missing track", 21, -0.5, 20.5);
123 fList->AddAtAndExpand(hNumberOfAdditionalTracks, kNumberOfAdditionalTracks);
124
125 // number of clusters
126 TH1F *hNumberOfClusters = new TH1F("hNumberOfClusters", "nb of clusters /track", 21, -0.5, 20.5);
127 fList->AddAtAndExpand(hNumberOfClusters, kNumberOfClusters);
128 TH1F *hNumberOfClustersM = new TH1F("hNumberOfClustersM", "nb of clusters /matched track", 21, -0.5, 20.5);
129 fList->AddAtAndExpand(hNumberOfClustersM, kNumberOfClustersM);
130 TH1F *hNumberOfClustersF = new TH1F("hNumberOfClustersF", "nb of clusters /fake track", 21, -0.5, 20.5);
131 fList->AddAtAndExpand(hNumberOfClustersF, kNumberOfClustersF);
132 TH1F *hNumberOfClustersMC = new TH1F("hNumberOfClustersMC", "nb of clusters /MC track", 21, -0.5, 20.5);
133 fList->AddAtAndExpand(hNumberOfClustersMC, kNumberOfClustersMC);
134 TH1F *hFractionOfMatchedClusters = new TH1F("hFractionOfMatchedClusters", "nb of matched clusters / nb of clusters", 110, 0., 1.1);
135 fList->AddAtAndExpand(hFractionOfMatchedClusters, kFractionOfMatchedClusters);
136 TH1F *hFractionOfConnectedClusters = new TH1F("hFractionOfConnectedClusters", "nb of connected clusters / nb of clusters in fake tracks", 110, 0., 1.1);
137 fList->AddAtAndExpand(hFractionOfConnectedClusters, kFractionOfConnectedClusters);
138
139 // number of fired chambers
140 TH1F *hNumberOfChamberHit = new TH1F("hNumberOfChamberHit", "nb of chambers hit /track", 16, -0.5, 15.5);
141 fList->AddAtAndExpand(hNumberOfChamberHit, kNumberOfChamberHit);
142 TH1F *hNumberOfChamberHitM = new TH1F("hNumberOfChamberHitM", "nb of chambers hit /matched track", 16, -0.5, 15.5);
143 fList->AddAtAndExpand(hNumberOfChamberHitM, kNumberOfChamberHitM);
144 TH1F *hNumberOfChamberHitF = new TH1F("hNumberOfChamberHitF", "nb of chambers hit /fake track", 16, -0.5, 15.5);
145 fList->AddAtAndExpand(hNumberOfChamberHitF, kNumberOfChamberHitF);
146
147 // chi2
148 TH1F *hChi2PerDof = new TH1F("hChi2PerDof", "track chi2/d.o.f.", 100, 0., 20.);
149 fList->AddAtAndExpand(hChi2PerDof, kChi2PerDof);
150 TH1F *hChi2PerDofM = new TH1F("hChi2PerDofM", "matched track chi2/d.o.f.", 100, 0., 20.);
151 fList->AddAtAndExpand(hChi2PerDofM, kChi2PerDofM);
152 TH1F *hChi2PerDofF = new TH1F("hChi2PerDofF", "fake track chi2/d.o.f.", 100, 0., 20.);
153 fList->AddAtAndExpand(hChi2PerDofF, kChi2PerDofF);
154
155 // chi2 versus number of clusters
156 TH2F *hChi2PerDofVsNClusters = new TH2F("hChi2PerDofVsNClusters", "track chi2/d.o.f. versus nb of clusters", 21, -0.5, 20.5, 100, 0., 20.);
157 fList->AddAtAndExpand(hChi2PerDofVsNClusters, kChi2PerDofVsNClusters);
158 TH2F *hChi2PerDofVsNClustersM = new TH2F("hChi2PerDofVsNClustersM", "matched track chi2/d.o.f. versus nb of clusters", 21, -0.5, 20.5, 100, 0., 20.);
159 fList->AddAtAndExpand(hChi2PerDofVsNClustersM, kChi2PerDofVsNClustersM);
160 TH2F *hChi2PerDofVsNClustersF = new TH2F("hChi2PerDofVsNClustersF", "fake track chi2/d.o.f. versus nb of clusters", 21, -0.5, 20.5, 100, 0., 20.);
161 fList->AddAtAndExpand(hChi2PerDofVsNClustersF, kChi2PerDofVsNClustersF);
162
163 // chi2 versus number of fired chambers
164 TH2F *hChi2PerDofVsNChamberHit = new TH2F("hChi2PerDofVsNChamberHit", "track chi2/d.o.f. versus nb of fired chambers", 16, -0.5, 15.5, 100, 0., 20.);
165 fList->AddAtAndExpand(hChi2PerDofVsNChamberHit, kChi2PerDofVsNChamberHit);
166 TH2F *hChi2PerDofVsNChamberHitM = new TH2F("hChi2PerDofVsNChamberHitM", "matched track chi2/d.o.f. versus nb of fired chambers", 16, -0.5, 15.5, 100, 0., 20.);
167 fList->AddAtAndExpand(hChi2PerDofVsNChamberHitM, kChi2PerDofVsNChamberHitM);
168 TH2F *hChi2PerDofVsNChamberHitF = new TH2F("hChi2PerDofVsNChamberHitF", "fake track chi2/d.o.f. versus nb of fired chambers", 16, -0.5, 15.5, 100, 0., 20.);
169 fList->AddAtAndExpand(hChi2PerDofVsNChamberHitF, kChi2PerDofVsNChamberHitF);
170
171 // physics quantities
172 TH1F *hP = new TH1F("hP", "Muon P distribution (GeV/c)", 100, 0., 200.);
173 fList->AddAtAndExpand(hP, kP);
174 TH1F *hPM = new TH1F("hPM", "matched track P distribution (GeV/c)", 100, 0., 200.);
175 fList->AddAtAndExpand(hPM, kPM);
176 TH1F *hPF = new TH1F("hPF", "fake track P distribution (GeV/c)", 100, 0., 200.);
177 fList->AddAtAndExpand(hPF, kPF);
178 TH1F *hPt = new TH1F("hPt", "Muon Pt distribution (GeV/c)", 100, 0., 20.);
179 fList->AddAtAndExpand(hPt, kPt);
180 TH1F *hPtM = new TH1F("hPtM", "matched track Pt distribution (GeV/c)", 100, 0., 20.);
181 fList->AddAtAndExpand(hPtM, kPtM);
182 TH1F *hPtF = new TH1F("hPtF", "fake track Pt distribution (GeV/c)", 100, 0., 20.);
183 fList->AddAtAndExpand(hPtF, kPtF);
184 TH1F *hEta = new TH1F("hEta", "Muon pseudo-rapidity distribution", 100, -10., 0.);
185 fList->AddAtAndExpand(hEta , kEta );
186 TH1F *hEtaM = new TH1F("hEtaM", "matched track pseudo-rapidity distribution", 100, -10., 0.);
187 fList->AddAtAndExpand(hEtaM, kEtaM);
188 TH1F *hEtaF = new TH1F("hEtaF", "fake track pseudo-rapidity distribution", 100, -10., 0.);
189 fList->AddAtAndExpand(hEtaF, kEtaF);
190 TH1F *hPhi = new TH1F("hPhi", "Muon phi distribution", 100, -1., 9.);
191 fList->AddAtAndExpand(hPhi, kPhi);
192 TH1F *hPhiM = new TH1F("hPhiM", "matched track phi distribution", 100, -1., 9.);
193 fList->AddAtAndExpand(hPhiM, kPhiM);
194 TH1F *hPhiF = new TH1F("hPhiF", "fake track phi distribution", 100, -1., 9.);
195 fList->AddAtAndExpand(hPhiF, kPhiF);
196 TH1F *hDCA = new TH1F("hDCA", "Muon DCA distribution", 125, 0., 500.);
197 fList->AddAtAndExpand(hDCA, kDCA);
198 TH1F *hDCAM = new TH1F("hDCAM", "matched track DCA distribution", 125, 0., 500.);
199 fList->AddAtAndExpand(hDCAM, kDCAM);
200 TH1F *hDCAF = new TH1F("hDCAF", "fake track DCA distribution", 125, 0., 500.);
201 fList->AddAtAndExpand(hDCAF, kDCAF);
202
203 // global counters of tracks:
204 // - reconstructible = number of reconstructible tracks
205 // - reconstructed = number of reconstructed tracks
206 // - matched = number of reconstructed tracks matched with a simulated one (reconstructible or not)
207 // - matchedyet = number of reconstructed tracks matched with a simulated one that is not reconstructible
208 // - fake = number of fake tracks
209 // - connected = number of fake tracks connected to a reconstructible simulated track
210 // - additional = number of additional (fake) tracks compared to the number of reconstructible ones
211 fTrackCounters = new AliCounterCollection("trackCounters");
212 fTrackCounters->AddRubric("track", "reconstructible/reconstructed/matched/matchedyet/fake/connected/additional");
213 fTrackCounters->AddRubric("run", 1000000);
214 fTrackCounters->AddRubric("trig", "yes/no/unknown");
215 fTrackCounters->AddRubric("selected", "yes/no");
216 fTrackCounters->AddRubric("acc", "in/out/unknown");
217 fTrackCounters->Init();
218
219 // detailled counters of fake tracks:
220 fFakeTrackCounters = new AliCounterCollection("fakeTrackCounters");
221 fFakeTrackCounters->AddRubric("track", "fake/connected/additional/matchedyet/fake?");
222 fFakeTrackCounters->AddRubric("run", 1000000);
223 fFakeTrackCounters->AddRubric("file", 1000000);
224 fFakeTrackCounters->AddRubric("event", 1000000);
225 fFakeTrackCounters->AddRubric("trig", "yes/no/unknown");
226 fFakeTrackCounters->AddRubric("selected", "yes/no");
227 fFakeTrackCounters->AddRubric("acc", "in/out/unknown");
228 fFakeTrackCounters->Init();
229
230 // counters of tracks matched by position or by using MC labels
231 fMatchedTrackCounters = new AliCounterCollection("matchedTrackCounters");
232 fMatchedTrackCounters->AddRubric("position", "match/not match");
233 fMatchedTrackCounters->AddRubric("label", "match/not match/match other");
234 fMatchedTrackCounters->AddRubric("run", 1000000);
235 fMatchedTrackCounters->AddRubric("trig", "yes/no");
236 fMatchedTrackCounters->AddRubric("selected", "yes/no");
237 fMatchedTrackCounters->AddRubric("acc", "in/out");
238 fMatchedTrackCounters->Init();
239
240 // global counters of events
241 // - any = total number of events with reconstructed tracks
242 // - fake = number of events with fake track(s)
243 // - notconnected = number of events with fake tracks that are not connected to a reconstructible simulated track
244 // - additional = number of events with additional (fake) tracks compared to the number of reconstructible ones
245 // - matchedyet = number of events with reconstructed tracks matched with a simulated one that is not reconstructible
246 // if trig = yes: only the tracks matched with the trigger are considered in the above logic
247 fEventCounters = new AliCounterCollection("eventCounters");
248 fEventCounters->AddRubric("event", "any/fake/notconnected/additional/matchedyet");
249 fEventCounters->AddRubric("run", 1000000);
250 fEventCounters->AddRubric("trig", "any/yes");
251 fEventCounters->AddRubric("selected", "yes/no");
252 fEventCounters->Init();
253
254 // Disable printout of AliMCEvent
255 AliLog::SetClassDebugLevel("AliMCEvent",-1);
256
257 // Post data at least once per task to ensure data synchronisation (required for merging)
258 PostData(1, fList);
259 PostData(2, fTrackCounters);
260 PostData(3, fFakeTrackCounters);
261 PostData(4, fMatchedTrackCounters);
262 PostData(5, fEventCounters);
263}
264
265//________________________________________________________________________
266void AliAnalysisTaskMuonFakes::UserExec(Option_t *)
267{
268 /// Process event: looks for fakes...
269
270 // check that reconstructions parameters for that run have been properly set
271 if (fSigmaCut < 0) return;
272
273 // check physics selection
274 TString selected = (fInputHandler && fInputHandler->IsEventSelected() != 0) ? "selected:yes" : "selected:no";
275
276 // current file name
277 fCurrentFileName = CurrentFileName();
278 fCurrentFileName.ReplaceAll("alien://","");
279 fCurrentFileName.ReplaceAll("/","\\");
280 fCurrentFileName.ReplaceAll(":",";");
281
282 // Load ESD event
283 AliESDEvent* esd = dynamic_cast<AliESDEvent*>(InputEvent());
284
285 // Load MC event
286 AliMCEventHandler *mcH = 0;
287 if(MCEvent()) mcH = static_cast<AliMCEventHandler*>((AliAnalysisManager::GetAnalysisManager())->GetMCtruthEventHandler());
288
289 // get reconstructed and simulated tracks
290 AliMUONRecoCheck rc(esd,mcH);
291 AliMUONVTrackStore* muonTrackStore = rc.ReconstructedTracks(-1, kFALSE);
292 AliMUONVTrackStore* trackRefStore = rc.TrackRefs(-1);
293 if (!muonTrackStore || !trackRefStore) return;
294
295 // count the number of reconstructible tracks
296 TIter next(trackRefStore->CreateIterator());
297 AliMUONTrack* trackRef;
298 while ( ( trackRef = static_cast<AliMUONTrack*>(next()) ) ) {
299 if (trackRef->IsValid(fRequestedStationMask, fRequest2ChInSameSt45))
300 fTrackCounters->Count(Form("track:reconstructible/run:%d/trig:unknown/%s/acc:unknown", fCurrentRunNumber, selected.Data()));
301 }
302
303 // loop over ESD tracks
304 Int_t nTrackerTracks = 0;
305 Bool_t containTrack[2] = {kFALSE, kFALSE};
306 Bool_t containFakeTrack[2] = {kFALSE, kFALSE};
307 Bool_t containMatchedYetTrack[2] = {kFALSE, kFALSE};
308 AliMUONVTrackStore *fakeTrackStore = AliMUONESDInterface::NewTrackStore();
309 Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ;
310 for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
311
312 AliESDMuonTrack* esdTrack = esd->GetMuonTrack(iTrack);
313
314 // skip ghosts
315 if (!esdTrack->ContainTrackerData()) continue;
316 nTrackerTracks++;
317 containTrack[0] = kTRUE;
318
319 // trigger condition
320 Bool_t trigger = esdTrack->ContainTriggerData();
321 TString trig = "trig:";
322 if (trigger) {
323 trig += "yes";
324 containTrack[1] = kTRUE;
325 } else trig += "no";
326
327 // acceptance condition
328 TString acc = "acc:";
329 Double_t thetaTrackAbsEnd = TMath::ATan(esdTrack->GetRAtAbsorberEnd()/505.) * TMath::RadToDeg();
330 Double_t eta = esdTrack->Eta();
331 if (thetaTrackAbsEnd >= 2. && thetaTrackAbsEnd <= 9. && eta >= -4. && eta <= -2.5) acc += "in";
332 else acc += "out";
333
334 // fill global counters
335 fTrackCounters->Count(Form("track:reconstructed/run:%d/%s/%s/%s", fCurrentRunNumber, trig.Data(), selected.Data(), acc.Data()));
336
337 // find the corresponding MUON track
338 AliMUONTrack* muonTrack = static_cast<AliMUONTrack*>(muonTrackStore->FindObject(esdTrack->GetUniqueID()));
339
340 // get track info
341 Int_t nClusters = esdTrack->GetNClusters();
342 Int_t nChamberHit = 0;
343 for (Int_t ich=0; ich<10; ich++) if (esdTrack->IsInMuonClusterMap(ich)) nChamberHit++;
344 Double_t normalizedChi2 = esdTrack->GetChi2() / (2. * esdTrack->GetNHit() - 5);
345 Double_t p = esdTrack->P();
346 Double_t pT = esdTrack->Pt();
347 Double_t phi = esdTrack->Phi();
348 Double_t dca = esdTrack->GetDCA();
349
350 // fill global histograms
351 ((TH1F*)fList->UncheckedAt(kNumberOfClusters))->Fill(nClusters);
352 ((TH1F*)fList->UncheckedAt(kNumberOfChamberHit))->Fill(nChamberHit);
353 ((TH1F*)fList->UncheckedAt(kChi2PerDof))->Fill(normalizedChi2);
354 ((TH1F*)fList->UncheckedAt(kP))->Fill(p);
355 ((TH1F*)fList->UncheckedAt(kPt))->Fill(pT);
356 ((TH1F*)fList->UncheckedAt(kEta))->Fill(eta);
357 ((TH1F*)fList->UncheckedAt(kPhi))->Fill(phi);
358 ((TH1F*)fList->UncheckedAt(kDCA))->Fill(dca);
359 ((TH1F*)fList->UncheckedAt(kChi2PerDofVsNClusters))->Fill(nClusters,normalizedChi2);
360 ((TH1F*)fList->UncheckedAt(kChi2PerDofVsNChamberHit))->Fill(nChamberHit,normalizedChi2);
361
362 // try to match, by position, the reconstructed track with a simulated one
363 Int_t nMatchClustersByPosition = 0;
364 AliMUONTrack* matchedTrackRefByPosition = rc.FindCompatibleTrack(*muonTrack, *trackRefStore, nMatchClustersByPosition, kFALSE, fSigmaCut);
365 Int_t MCLabelByPosition = (matchedTrackRefByPosition) ? matchedTrackRefByPosition->GetUniqueID() : -1;
366
367 // try to match, by using MC labels, the reconstructed track with a simulated one
368 Int_t nMatchClustersByLabel = 0;
369 AliMUONTrack* matchedTrackRefByLabel = rc.FindCompatibleTrack(*muonTrack, *trackRefStore, nMatchClustersByLabel, kTRUE, fSigmaCut);
370 Int_t MCLabelByLabel = (matchedTrackRefByLabel) ? matchedTrackRefByLabel->GetUniqueID() : -1;
371
372 // fill global counters
373 TString positionCase = (MCLabelByPosition >= 0) ? "position:match" : "position:not match";
374 TString labelCase = "label:";
375 if (MCLabelByLabel >= 0 && MCLabelByPosition >= 0 && MCLabelByLabel != MCLabelByPosition) labelCase += "match other";
376 else if (MCLabelByLabel >= 0) labelCase += "match";
377 else labelCase += "not match";
378 fMatchedTrackCounters->Count(Form("%s/%s/run:%d/%s/%s/%s", positionCase.Data(), labelCase.Data(), fCurrentRunNumber, trig.Data(), selected.Data(), acc.Data()));
379 if ((MCLabelByLabel >= 0 && MCLabelByPosition < 0) || (MCLabelByLabel < 0 && MCLabelByPosition >= 0))
380 fFakeTrackCounters->Count(Form("track:fake?/run:%d/file:%s/event:%d/%s/%s/%s", fCurrentRunNumber, fCurrentFileName.Data(),
381 esd->GetEventNumberInFile(), trig.Data(), selected.Data(), acc.Data()));
382
383 // take actions according to the matching result we are interested in
384 Int_t nMatchClusters = (fUseLabel) ? nMatchClustersByLabel : nMatchClustersByPosition;
385 AliMUONTrack* matchedTrackRef = (fUseLabel) ? matchedTrackRefByLabel : matchedTrackRefByPosition;
386 if (matchedTrackRef) {
387
388 // fill global counters
389 fTrackCounters->Count(Form("track:matched/run:%d/%s/%s/%s", fCurrentRunNumber, trig.Data(), selected.Data(), acc.Data()));
390
391 // track matched with a trackRef that is not reconstructible
392 if (!matchedTrackRef->IsValid(fRequestedStationMask, fRequest2ChInSameSt45)) {
393
394 containMatchedYetTrack[0] = kTRUE;
395 if (trigger) containMatchedYetTrack[1] = kTRUE;
396
397 // fill global counters
398 fTrackCounters->Count(Form("track:matchedyet/run:%d/%s/%s/%s", fCurrentRunNumber, trig.Data(), selected.Data(), acc.Data()));
399 fFakeTrackCounters->Count(Form("track:matchedyet/run:%d/file:%s/event:%d/%s/%s/%s", fCurrentRunNumber, fCurrentFileName.Data(),
400 esd->GetEventNumberInFile(), trig.Data(), selected.Data(), acc.Data()));
401 }
402
403 // fill histograms
404 ((TH1F*)fList->UncheckedAt(kFractionOfMatchedClusters))->Fill(((Float_t) nMatchClusters) / ((Float_t) nClusters));
405 ((TH1F*)fList->UncheckedAt(kNumberOfClustersMC))->Fill(matchedTrackRef->GetNClusters());
406 ((TH1F*)fList->UncheckedAt(kNumberOfClustersM))->Fill(nClusters);
407 ((TH1F*)fList->UncheckedAt(kNumberOfChamberHitM))->Fill(nChamberHit);
408 ((TH1F*)fList->UncheckedAt(kChi2PerDofM))->Fill(normalizedChi2);
409 ((TH1F*)fList->UncheckedAt(kPM))->Fill(p);
410 ((TH1F*)fList->UncheckedAt(kPtM))->Fill(pT);
411 ((TH1F*)fList->UncheckedAt(kEtaM))->Fill(eta);
412 ((TH1F*)fList->UncheckedAt(kPhiM))->Fill(phi);
413 ((TH1F*)fList->UncheckedAt(kDCAM))->Fill(dca);
414 ((TH1F*)fList->UncheckedAt(kChi2PerDofVsNClustersM))->Fill(nClusters,normalizedChi2);
415 ((TH1F*)fList->UncheckedAt(kChi2PerDofVsNChamberHitM))->Fill(nChamberHit,normalizedChi2);
416
417 // remove already matched trackRefs
418 trackRefStore->Remove(*matchedTrackRef);
419
420 } else {
421
422 containFakeTrack[0] = kTRUE;
423 if (trigger) containFakeTrack[1] = kTRUE;
424
425 // fill global counters
426 fTrackCounters->Count(Form("track:fake/run:%d/%s/%s/%s", fCurrentRunNumber, trig.Data(), selected.Data(), acc.Data()));
427 fFakeTrackCounters->Count(Form("track:fake/run:%d/file:%s/event:%d/%s/%s/%s", fCurrentRunNumber, fCurrentFileName.Data(),
428 esd->GetEventNumberInFile(), trig.Data(), selected.Data(), acc.Data()));
429
430 // fill histograms
431 ((TH1F*)fList->UncheckedAt(kNumberOfClustersF))->Fill(nClusters);
432 ((TH1F*)fList->UncheckedAt(kNumberOfChamberHitF))->Fill(nChamberHit);
433 ((TH1F*)fList->UncheckedAt(kChi2PerDofF))->Fill(normalizedChi2);
434 ((TH1F*)fList->UncheckedAt(kPF))->Fill(p);
435 ((TH1F*)fList->UncheckedAt(kPtF))->Fill(pT);
436 ((TH1F*)fList->UncheckedAt(kEtaF))->Fill(eta);
437 ((TH1F*)fList->UncheckedAt(kPhiF))->Fill(phi);
438 ((TH1F*)fList->UncheckedAt(kDCAF))->Fill(dca);
439 ((TH1F*)fList->UncheckedAt(kChi2PerDofVsNClustersF))->Fill(nClusters,normalizedChi2);
440 ((TH1F*)fList->UncheckedAt(kChi2PerDofVsNChamberHitF))->Fill(nChamberHit,normalizedChi2);
441
442 // store fake tracks
443 fakeTrackStore->Add(*muonTrack);
444
445 }
446
447 } // end of loop over ESD tracks
448
449 // fill histogram and global counters
450 ((TH1F*)fList->UncheckedAt(kNumberOfTracks))->Fill(nTrackerTracks);
451 if (containTrack[0]) fEventCounters->Count(Form("event:any/run:%d/trig:any/%s", fCurrentRunNumber, selected.Data()));
452 if (containTrack[1]) fEventCounters->Count(Form("event:any/run:%d/trig:yes/%s", fCurrentRunNumber, selected.Data()));
453 if (containFakeTrack[0]) fEventCounters->Count(Form("event:fake/run:%d/trig:any/%s", fCurrentRunNumber, selected.Data()));
454 if (containFakeTrack[1]) fEventCounters->Count(Form("event:fake/run:%d/trig:yes/%s", fCurrentRunNumber, selected.Data()));
455 if (containMatchedYetTrack[0]) fEventCounters->Count(Form("event:matchedyet/run:%d/trig:any/%s", fCurrentRunNumber, selected.Data()));
456 if (containMatchedYetTrack[1]) fEventCounters->Count(Form("event:matchedyet/run:%d/trig:yes/%s", fCurrentRunNumber, selected.Data()));
457
458 // count the number of not connected and additional fake tracks
459 if (fakeTrackStore->GetSize() > 0) {
460
461 // remove the most connected fake tracks
462 Int_t nFreeMissingTracks = RemoveConnectedFakes(*fakeTrackStore, *trackRefStore);
463
464 if (fakeTrackStore->GetSize() > 0) {
465
466 // fill global counters
467 fEventCounters->Count(Form("event:notconnected/run:%d/trig:any/%s", fCurrentRunNumber, selected.Data()));
468
469 // check status of remaining fakes with respect to the matching with trigger
470 Bool_t containMatchedFake = kFALSE;
471 Bool_t containUnmatchedFake = kFALSE;
472 AliMUONTrack* fakeTrack = 0x0;
473 TIter next3(fakeTrackStore->CreateIterator());
474 while ( ( fakeTrack = static_cast<AliMUONTrack*>(next3()) ) ) {
475 if (fakeTrack->GetMatchTrigger() > 0) containMatchedFake = kTRUE;
476 else containUnmatchedFake = kTRUE;
477 }
478
479 // fill global counters
480 if (containMatchedFake) fEventCounters->Count(Form("event:notconnected/run:%d/trig:yes/%s", fCurrentRunNumber, selected.Data()));
481
482 // remove the remaining free reconstructible tracks
483 Int_t nAdditionalTracks = fakeTrackStore->GetSize() - nFreeMissingTracks;
484
485 if (nAdditionalTracks > 0) {
486
487 // fill histogram and global counters
488 ((TH1F*)fList->UncheckedAt(kNumberOfAdditionalTracks))->Fill(nAdditionalTracks);
489 fEventCounters->Count(Form("event:additional/run:%d/trig:any/%s", fCurrentRunNumber, selected.Data()));
490 if (!containUnmatchedFake) { // all matched
491 fTrackCounters->Count(Form("track:additional/run:%d/trig:yes/%s/acc:unknown", fCurrentRunNumber, selected.Data()), nAdditionalTracks);
492 fFakeTrackCounters->Count(Form("track:additional/run:%d/file:%s/event:%d/trig:yes/%s/acc:unknown", fCurrentRunNumber, fCurrentFileName.Data(), esd->GetEventNumberInFile(), selected.Data()), nAdditionalTracks);
493 fEventCounters->Count(Form("event:additional/run:%d/trig:yes/%s", fCurrentRunNumber, selected.Data()));
494 } else if (!containMatchedFake) { // none matched
495 fTrackCounters->Count(Form("track:additional/run:%d/trig:no/%s/acc:unknown", fCurrentRunNumber, selected.Data()), nAdditionalTracks);
496 fFakeTrackCounters->Count(Form("track:additional/run:%d/file:%s/event:%d/trig:no/%s/acc:unknown", fCurrentRunNumber, fCurrentFileName.Data(), esd->GetEventNumberInFile(), selected.Data()), nAdditionalTracks);
497 } else { // mixed
498 fTrackCounters->Count(Form("track:additional/run:%d/trig:unknown/%s/acc:unknown", fCurrentRunNumber, selected.Data()), nAdditionalTracks);
499 fFakeTrackCounters->Count(Form("track:additional/run:%d/file:%s/event:%d/trig:unknown/%s/acc:unknown", fCurrentRunNumber, fCurrentFileName.Data(), esd->GetEventNumberInFile(), selected.Data()), nAdditionalTracks);
500 fEventCounters->Count(Form("event:additional/run:%d/trig:yes/%s", fCurrentRunNumber, selected.Data()));
501 }
502
503 }
504
505 }
506
507 }
508
509 delete fakeTrackStore;
510
511 // Post final data
512 PostData(1, fList);
513 PostData(2, fTrackCounters);
514 PostData(3, fFakeTrackCounters);
515 PostData(4, fMatchedTrackCounters);
516 PostData(5, fEventCounters);
517}
518
519//________________________________________________________________________
520void AliAnalysisTaskMuonFakes::NotifyRun()
521{
522 /// Prepare processing of new run: load corresponding OCDB objects...
523
524 // load necessary data from OCDB
525 AliCDBManager::Instance()->SetDefaultStorage(fRecoParamLocation.Data());
526 AliCDBManager::Instance()->SetRun(fCurrentRunNumber);
527 AliMUONRecoParam* recoParam = AliMUONCDB::LoadRecoParam();
528 if (!recoParam) {
529 fRequestedStationMask = 0;
530 fRequest2ChInSameSt45 = kFALSE;
531 fSigmaCut = -1.;
532 AliError("--> skip this run");
533 return;
534 }
535
536 // compute the mask of requested stations from recoParam
537 fRequestedStationMask = 0;
538 for (Int_t i = 0; i < 5; i++) if (recoParam->RequestStation(i)) fRequestedStationMask |= ( 1 << i );
539
540 // get from recoParam whether a track need 2 chambers hit in the same station (4 or 5) or not to be reconstructible
541 fRequest2ChInSameSt45 = !recoParam->MakeMoreTrackCandidates();
542
543 // get sigma cut from recoParam to associate clusters with TrackRefs in case the labels are not used
544 fSigmaCut = (recoParam->ImproveTracks()) ? recoParam->GetSigmaCutForImprovement() : recoParam->GetSigmaCutForTracking();
545}
546
547//________________________________________________________________________
548void AliAnalysisTaskMuonFakes::Terminate(Option_t *)
549{
550 /// Draw results to the screen and print statistics.
551
552 // recover output objects
553 fList = static_cast<TObjArray*> (GetOutputData(1));
554 if (!fList) return;
555 fTrackCounters = static_cast<AliCounterCollection*> (GetOutputData(2));
556 fFakeTrackCounters = static_cast<AliCounterCollection*> (GetOutputData(3));
557 fMatchedTrackCounters = static_cast<AliCounterCollection*> (GetOutputData(4));
558 fEventCounters = static_cast<AliCounterCollection*> (GetOutputData(5));
559
560 // add canvas to compare histograms
561 fCanvases = new TObjArray(1000);
562 fCanvases->SetOwner();
563 TCanvas *cFakesSummary1 = new TCanvas("cFakesSummary1","cFakesSummary1",1200,600);
564 fCanvases->AddAtAndExpand(cFakesSummary1, 0);
565 TCanvas *cFakesSummary2 = new TCanvas("cFakesSummary2","cFakesSummary2",1200,600);
566 fCanvases->AddAtAndExpand(cFakesSummary2, 1);
567
568 // display
569 Int_t iHist1[8] = {kNumberOfClusters, kChi2PerDof, kP, kEta, kNumberOfChamberHit, kDCA, kPt, kPhi};
570 cFakesSummary1->Divide(4,2);
571 for (Int_t i=0; i<8; i++) {
572 cFakesSummary1->cd(i+1);
573 cFakesSummary1->GetPad(i+1)->SetLogy();
574 ((TH1F*)fList->UncheckedAt(iHist1[i]))->SetMinimum(0.5);
575 ((TH1F*)fList->UncheckedAt(iHist1[i]))->DrawCopy();
576 ((TH1F*)fList->UncheckedAt(iHist1[i]+1))->SetLineColor(4);
577 ((TH1F*)fList->UncheckedAt(iHist1[i]+1))->DrawCopy("sames");
578 ((TH1F*)fList->UncheckedAt(iHist1[i]+2))->SetLineColor(2);
579 ((TH1F*)fList->UncheckedAt(iHist1[i]+2))->SetFillColor(2);
580 ((TH1F*)fList->UncheckedAt(iHist1[i]+2))->SetFillStyle(3017);
581 ((TH1F*)fList->UncheckedAt(iHist1[i]+2))->DrawCopy("sames");
582 }
583
584 Int_t iHist2[2] = {kChi2PerDofVsNClusters, kChi2PerDofVsNChamberHit};
585 cFakesSummary2->Divide(2);
586 for (Int_t i=0; i<2; i++) {
587 cFakesSummary2->cd(i+1);
588 ((TH2F*)fList->UncheckedAt(iHist2[i]+1))->SetMarkerColor(4);
589 ((TH2F*)fList->UncheckedAt(iHist2[i]+1))->DrawCopy();
590 ((TH2F*)fList->UncheckedAt(iHist2[i]+2))->SetMarkerColor(2);
591 ((TH2F*)fList->UncheckedAt(iHist2[i]+2))->SetMarkerStyle(7);
592 ((TH2F*)fList->UncheckedAt(iHist2[i]+2))->DrawCopy("sames");
593 }
594
595 // print
596 if (fTrackCounters && fFakeTrackCounters && fMatchedTrackCounters && fEventCounters) {
597 printf("\nGlobal statistics of reconstructed tracks matched or not with the trigger:\n");
598 fTrackCounters->Print("track/trig");
599 printf("\nGlobal statistics of pathological tracks matched or not with the trigger:\n");
600 fFakeTrackCounters->Print("track/trig");
601 printf("\nDetailled statistics of tracks matched per label vs position:\n");
602 fMatchedTrackCounters->Print("label/position");
603 printf("\nGlobal statistics of events containing pathological tracks:\n");
604 fEventCounters->Print("event/trig");
605 }
606
607 printf("\nREMINDER: results are relevent provided that you use the same recoParams as for the reconstruction\n");
608}
609
610//________________________________________________________________________
611Int_t AliAnalysisTaskMuonFakes::RemoveConnectedFakes(AliMUONVTrackStore &fakeTrackStore, AliMUONVTrackStore &trackRefStore)
612{
613 /// loop over reconstructible TrackRef not associated with reconstructed track:
614 /// for each of them, find and remove the most connected the fake track, if any,
615 /// and fill the histograms with the fraction of connected clusters.
616 /// Return the number of reconstructible track not connected to any fake
617
618 Int_t nFreeMissingTracks = 0;
619
620 // loop over trackRefs
621 TIter next(trackRefStore.CreateIterator());
622 AliMUONTrack* trackRef;
623 while ( ( trackRef = static_cast<AliMUONTrack*>(next()) ) ) {
624
625 // skip not reconstructible trackRefs
626 if (!trackRef->IsValid(fRequestedStationMask, fRequest2ChInSameSt45)) continue;
627
628 Int_t label = trackRef->GetUniqueID();
629
630 // look for the most connected fake track
631 AliMUONTrack *connectedFake = 0x0;
632 Double_t fractionOfConnectedClusters = 0.;
633 TIter next2(fakeTrackStore.CreateIterator());
634 AliMUONTrack* fakeTrack;
635 while ( ( fakeTrack = static_cast<AliMUONTrack*>(next2()) ) ) {
636
637 // get the number of connected clusters
638 Int_t nConnectedClusters = 0;
639 if (fUseLabel) { // by using the MC label
640 for (Int_t iCl = 0; iCl < fakeTrack->GetNClusters(); iCl++)
641 if (((AliMUONTrackParam*) fakeTrack->GetTrackParamAtCluster()->UncheckedAt(iCl))->GetClusterPtr()->GetMCLabel() == label)
642 nConnectedClusters++;
643 } else { // by comparing cluster/TrackRef positions
644 Bool_t compTrack[10];
645 nConnectedClusters = fakeTrack->FindCompatibleClusters(*trackRef, fSigmaCut, compTrack);
646 }
647
648 // skip non-connected fake tracks
649 if (nConnectedClusters == 0) continue;
650
651 // check if it is the most connected fake track
652 Double_t f = ((Double_t)nConnectedClusters) / ((Double_t)fakeTrack->GetNClusters());
653 if (f > fractionOfConnectedClusters) {
654 connectedFake = fakeTrack;
655 fractionOfConnectedClusters = f;
656 }
657
658 }
659
660 if (connectedFake) {
661
662 // find the corresponding ESD MUON track
663 AliESDEvent* esd = dynamic_cast<AliESDEvent*>(InputEvent());
664 TIter next3(static_cast<TClonesArray*>(esd->FindListObject("MuonTracks")));
665 AliESDMuonTrack* esdTrack = 0x0;
666 while ((esdTrack = static_cast<AliESDMuonTrack*>(next3())) && esdTrack->GetUniqueID() != connectedFake->GetUniqueID()) {}
667 if (!esdTrack) {
668 AliError("unable to find the corresponding ESD track???");
669 continue;
670 }
671
672 // trigger condition
673 TString trig = (esdTrack->ContainTriggerData()) ? "trig:yes" : "trig:no";
674
675 // acceptance condition
676 Double_t thetaTrackAbsEnd = TMath::ATan(esdTrack->GetRAtAbsorberEnd()/505.) * TMath::RadToDeg();
677 Double_t eta = esdTrack->Eta();
678 TString acc = (thetaTrackAbsEnd >= 2. && thetaTrackAbsEnd <= 9. && eta >= -4. && eta <= -2.5) ? "acc:in" : "acc:out";
679
680 // check physics selection
681 TString selected = (fInputHandler && fInputHandler->IsEventSelected()) ? "selected:yes" : "selected:no";
682
683 // fill histogram and counters
684 ((TH1F*)fList->UncheckedAt(kFractionOfConnectedClusters))->Fill(fractionOfConnectedClusters);
685 fTrackCounters->Count(Form("track:connected/run:%d/%s/%s/%s", fCurrentRunNumber, trig.Data(), selected.Data(), acc.Data()));
686 fFakeTrackCounters->Count(Form("track:connected/run:%d/file:%s/event:%d/%s/%s/%s", fCurrentRunNumber, fCurrentFileName.Data(),
687 esd->GetEventNumberInFile(), trig.Data(), selected.Data(), acc.Data()));
688
689 // remove the most connected fake track
690 fakeTrackStore.Remove(*connectedFake);
691
692 } else nFreeMissingTracks++;
693
694 }
695
696 return nFreeMissingTracks;
697
698}
699