]>
Commit | Line | Data |
---|---|---|
b220e595 | 1 | #include "TFile.h" |
2 | #include "TTree.h" | |
3 | #include "TChain.h" | |
4 | #include "TH1F.h" | |
5 | #include "TH2F.h" | |
6 | #include "TProfile.h" | |
7 | #include "TCanvas.h" | |
8 | ||
9 | #include "AliAnalysisManager.h" | |
10 | #include "AliESDEvent.h" | |
11 | #include "AliAODEvent.h" | |
12 | #include "AliMCEvent.h" | |
13 | #include "AliESDInputHandler.h" | |
14 | #include "AliAODHandler.h" | |
15 | #include "AliMCEventHandler.h" | |
16 | #include "AliLog.h" | |
17 | #include "AliESDTrdTrack.h" | |
18 | ||
19 | #include "AliTRDtrackletMCM.h" | |
20 | #include "AliVParticle.h" | |
21 | #include "AliMCParticle.h" | |
22 | #include "AliTrackReference.h" | |
23 | #include "AliMagF.h" | |
24 | #include "TGeoGlobalMagField.h" | |
25 | #include "AliTRDtrackletWord.h" | |
26 | ||
27 | #include "AliTRDonlineTrackletQA.h" | |
28 | ||
29 | ClassImp(AliTRDonlineTrackletQA) | |
30 | ||
31 | AliTRDonlineTrackletQA::AliTRDonlineTrackletQA(const char *name) : | |
32 | AliAnalysisTask(name, ""), | |
33 | fESD(0x0), | |
34 | fInputHandler(0x0), | |
35 | fInputEvent(0x0), | |
36 | fOutputAOD(0x0), | |
37 | fMCEvent(0x0), | |
38 | fTrackletsRaw(0x0), | |
39 | fTrackletsSim(0x0), | |
40 | fOutputList(0x0), | |
41 | fHistYpos(0x0), | |
42 | fHistYres(0x0), | |
43 | fHistYresDy(0x0), | |
44 | fHistdY(0x0), | |
45 | fHistdYres(0x0), | |
46 | fHistYresESD(0x0), | |
47 | fHistdYresESD(0x0), | |
48 | fHistCanddY(0x0), | |
49 | fHistFounddY(0x0), | |
50 | fHistTrklPerRef(0x0), | |
51 | fHistdYdYref(0x0), | |
52 | fHistYposRaw(0x0), | |
53 | fHistdYRaw(0x0), | |
54 | fHistYdYRaw(0x0), | |
55 | fHistZrow(0x0), | |
56 | fHistZrowRaw(0x0), | |
57 | fHistPid(0x0), | |
58 | fHistPidRaw(0x0), | |
59 | fHistYdiff(0x0), | |
60 | fHistdYdiff(0x0), | |
61 | fHistdYdYraw(0x0), | |
62 | fTreeTracklets(0x0), | |
63 | fMinPt(1.), | |
64 | fGeo(new AliTRDgeometry), | |
65 | fNevent(0), | |
66 | fTrackletFile(0x0), | |
67 | fTrackletTree(0x0), | |
68 | fTrackletTreeRaw(0x0) | |
69 | { | |
70 | // ctor | |
71 | ||
72 | DefineInput(0, TChain::Class()); | |
73 | DefineInput(1, TTree::Class()); | |
74 | ||
75 | DefineOutput(0, TTree::Class()); | |
76 | DefineOutput(1, TList::Class()); | |
77 | } | |
78 | ||
79 | AliTRDonlineTrackletQA::~AliTRDonlineTrackletQA() | |
80 | { | |
81 | // dtor | |
82 | ||
83 | delete fGeo; | |
84 | } | |
85 | ||
86 | void AliTRDonlineTrackletQA::ConnectInputData(const Option_t *option) | |
87 | { | |
88 | fInputHandler = (AliInputEventHandler*) AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler(); | |
89 | if (fInputHandler) | |
90 | fInputEvent = fInputHandler->GetEvent(); | |
91 | ||
92 | AliMCEventHandler *mcH = (AliMCEventHandler*) AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler(); | |
93 | if (mcH) | |
94 | fMCEvent = mcH->MCEvent(); | |
95 | } | |
96 | ||
97 | void AliTRDonlineTrackletQA::CreateOutputObjects() | |
98 | { | |
99 | OpenFile(1); | |
100 | ||
101 | fOutputList = new TList(); | |
102 | ||
103 | fHistYpos = new TH1F("ypos", | |
104 | "Tracklet (sim) y-position;y (cm);count", | |
105 | 8192/32, -4096*160e-4, 4095*160e-4); | |
106 | fHistYposRaw = new TH1F("ypos-raw", | |
107 | "Tracklet (raw) y-position;y (cm);count", | |
108 | 130, -65, 65); | |
109 | fHistYres = new TH1F("yres", | |
110 | "Tracklet (sim) #Deltay;y_{tracklet}-y_{MC} (cm);count", | |
111 | 8192/32, -4096/32*160e-4, 4095/32*160e-4); | |
112 | fHistYresDy = new TH2F("yresdy", | |
113 | "Tracklet (sim) #Deltay;y_{tracklet}-y_{MC} (cm);deflection (bin)", | |
114 | 8192/32, -4096/32*160e-4, 4095/32*160e-4, | |
115 | 128, -64.5, 63.5); | |
116 | fHistYresESD = new TH1F("yresesd", | |
117 | "Tracklet #Deltay;y (cm);count", | |
118 | 100, -10, 10); | |
119 | fHistYdiff = new TH1F("ydiff", | |
120 | "Tracklet #Deltay (sim - raw);y_{sim}-y_{raw} (cm);count", | |
121 | 8192/32, -4096/32*160e-4, 4095/32*160e-4); | |
122 | for (Int_t iLayer = 0; iLayer < 6; iLayer++) { | |
123 | fHistYlocal[iLayer] = new TH2F(Form("ylocal_%i", iLayer), | |
124 | Form("Tracklet local y, layer %i;y_{MC} (pad width);y_{trkl} (pad width)", iLayer), | |
125 | 100, -1, 1, 100, -1, 1); | |
126 | } | |
127 | ||
128 | fHistdY = new TH1F("dy", | |
129 | "deflection (sim);dy (140 #mum)", | |
130 | 128, -64.5, 63.5); | |
131 | fHistdYRaw = new TH1F("dy-raw", | |
132 | "deflection (sim);dy (140 #mum)", | |
133 | 128, -64.5, 63.5); | |
134 | fHistdYres = new TH1F("dyres", | |
135 | "deflection residual;dy (cm)", | |
136 | 128, -1., 1.); | |
137 | fHistdYresESD = new TH1F("dyresesd", | |
138 | "deflection residual;dy (cm)", | |
139 | 128, -1., 1.); | |
140 | fHistCanddY = new TH1F("dycand", | |
141 | "deflection;dy (140 #mum)", | |
142 | 128, -64.5, 63.5); | |
143 | fHistFounddY = new TH1F("dyfound", | |
144 | "deflection;dy (140 #mum)", | |
145 | 128, -64.5, 63.5); | |
146 | fHistdYdiff = new TH1F("dydiff", | |
147 | "deflection #Deltady;dy_{sim}-dy_{raw} (140 #mum)", | |
148 | 128, -64.5, 63.5); | |
149 | fHistdYdYraw = new TH2F("dydyraw", | |
150 | "deflection from sim. vs raw;dy_{sim} (140 #mum);dy_{raw} (140 #mum)", | |
151 | 128, -64.5, 63.5, 128, -64.5, 63.5); | |
152 | ||
153 | fHistTrklPerRef = new TH1F("trklperref", | |
154 | "No. of tracklets per track reference;no. of tracklets", | |
155 | 10, -0.5, 9.5); | |
156 | ||
157 | fHistdYdYref = new TH2F("dydyref", | |
158 | "deflection vs. deflection from track reference;dy_{ref} (140 #mum);dy (140 #mum)", | |
159 | 128, -64.5, 63.5, 128, -64.5, 63.5); | |
160 | ||
161 | fHistZrow = new TH1F("zrow", | |
162 | "z-position;pad row", | |
163 | 16, -0.5, 15.5); | |
164 | fHistZrowRaw = new TH1F("zrow-raw", | |
165 | "z-position;pad row", | |
166 | 16, -0.5, 15.5); | |
167 | ||
168 | fHistPid = new TH1F("pid", | |
169 | "pid", | |
170 | 256, -0.5, 255.5); | |
171 | fHistPidRaw = new TH1F("pid-raw", | |
172 | "pid", | |
173 | 256, -0.5, 255.5); | |
174 | ||
175 | fHistYdYRaw = new TH2F("ydyraw", | |
176 | "y vs dy (raw tracklets);y (cm);dy (140 #mum)", | |
177 | 8192/32, -4096*160e-4, 4095*160e-4, | |
178 | 128, -64.5, 63.5); | |
179 | ||
180 | fTreeTracklets = new TTree("trkl", "trkl"); | |
181 | fTreeTracklets->Branch("y", &fY); | |
182 | fTreeTracklets->Branch("dy", &fDY); | |
183 | fTreeTracklets->Branch("ydiff", &fYdiff); | |
184 | fTreeTracklets->Branch("dydiff", &fDYdiff); | |
185 | fTreeTracklets->Branch("q0", &fQ0); | |
186 | fTreeTracklets->Branch("q1", &fQ1); | |
187 | fTreeTracklets->Branch("nhits", &fNHits); | |
188 | ||
189 | fOutputList->Add(fHistYpos); | |
190 | fOutputList->Add(fHistdY); | |
191 | fOutputList->Add(fHistZrow); | |
192 | fOutputList->Add(fHistPid); | |
193 | ||
194 | fOutputList->Add(fHistYres); | |
195 | fOutputList->Add(fHistYresDy); | |
196 | fOutputList->Add(fHistCanddY); | |
197 | fOutputList->Add(fHistFounddY); | |
198 | fOutputList->Add(fHistTrklPerRef); | |
199 | fOutputList->Add(fHistdYres); | |
200 | fOutputList->Add(fHistYresESD); | |
201 | fOutputList->Add(fHistdYresESD); | |
202 | fOutputList->Add(fHistdYdYref); | |
203 | ||
204 | for (Int_t iLayer = 0; iLayer < 6; iLayer++) | |
205 | fOutputList->Add(fHistYlocal[iLayer]); | |
206 | ||
207 | fOutputList->Add(fHistYposRaw); | |
208 | fOutputList->Add(fHistdYRaw); | |
209 | fOutputList->Add(fHistZrowRaw); | |
210 | fOutputList->Add(fHistPidRaw); | |
211 | fOutputList->Add(fHistYdYRaw); | |
212 | ||
213 | fOutputList->Add(fHistYdiff); | |
214 | fOutputList->Add(fHistdYdiff); | |
215 | fOutputList->Add(fHistdYdYraw); | |
216 | ||
217 | fOutputList->Add(fTreeTracklets); | |
218 | } | |
219 | ||
220 | void AliTRDonlineTrackletQA::Exec(const Option_t *option) | |
221 | { | |
222 | // execute this for each event | |
223 | ||
224 | // connect input data | |
225 | fTrackletTree = (TTree*) GetInputData(1); | |
226 | if (!fTrackletTree) | |
227 | return; | |
228 | ||
229 | fTrackletTree->SetBranchAddress("tracklets_sim", &fTrackletsSim); | |
230 | fTrackletTree->SetBranchAddress("tracklets_raw", &fTrackletsRaw); | |
231 | fTrackletTree->GetEntry(fTrackletTree->GetEntriesFast()-1); | |
232 | ||
233 | // prepare raw tracklets for comparison | |
234 | Int_t detRaw; | |
235 | Int_t robRaw; | |
236 | Int_t mcmRaw; | |
237 | Int_t yRaw; | |
238 | Int_t dyRaw; | |
239 | TTree trklRaw("raw tracklets", "raw tracklets"); | |
240 | trklRaw.Branch("det", &detRaw); | |
241 | trklRaw.Branch("rob", &robRaw); | |
242 | trklRaw.Branch("mcm", &mcmRaw); | |
243 | trklRaw.Branch("y", &yRaw); | |
244 | trklRaw.Branch("dy", &dyRaw); | |
245 | trklRaw.SetDirectory(0x0); | |
246 | // prepare simulated tracklets for comparison | |
247 | Int_t detSim; | |
248 | Int_t robSim; | |
249 | Int_t mcmSim; | |
250 | Int_t ySim; | |
251 | Int_t dySim; | |
252 | TTree trklSim("sim tracklets", "sim tracklets"); | |
253 | trklSim.Branch("det", &detSim); | |
254 | trklSim.Branch("rob", &robSim); | |
255 | trklSim.Branch("mcm", &mcmSim); | |
256 | trklSim.Branch("y", &ySim); | |
257 | trklSim.Branch("dy", &dySim); | |
258 | trklSim.SetDirectory(0x0); | |
259 | ||
260 | // ----- simulated tracklets ----- | |
261 | AliTRDtrackletMCM *trkl = 0x0; | |
262 | if (fTrackletsSim) { | |
263 | for (Int_t iTracklet = 0; iTracklet < fTrackletsSim->GetEntries(); iTracklet++) { | |
264 | trkl = (AliTRDtrackletMCM*) (*fTrackletsSim)[iTracklet]; | |
265 | Int_t label = trkl->GetLabel(); | |
266 | // if (label > -1 && label < maxTracks) | |
267 | // mcTrackToTrackletMCM[label].idx[mcTrackToTrackletMCM[label].size < 10 ? mcTrackToTrackletMCM[label].size++ : 0] = iTracklet; | |
268 | fHistYpos->Fill(trkl->GetY()); | |
269 | fHistdY->Fill(trkl->GetdY()); | |
270 | fHistZrow->Fill(trkl->GetZbin()); | |
271 | fHistPid->Fill(trkl->GetPID()); | |
272 | ||
273 | detSim = trkl->GetDetector(); | |
274 | robSim = trkl->GetROB(); | |
275 | mcmSim = trkl->GetMCM(); | |
276 | ySim = trkl->GetYbin(); | |
277 | dySim = trkl->GetdY(); | |
278 | trklSim.Fill(); | |
279 | ||
280 | PlotMC(trkl); | |
281 | PlotESD(trkl); | |
282 | } | |
283 | } | |
284 | ||
285 | // ----- raw tracklets ----- | |
286 | if (fTrackletsRaw) { | |
287 | for (Int_t iTracklet = 0; iTracklet < fTrackletsRaw->GetEntries(); iTracklet++) { | |
288 | AliTRDtrackletWord *trklWord = (AliTRDtrackletWord*) (*fTrackletsRaw)[iTracklet]; | |
289 | // remove unwanted chambers | |
290 | if (trklWord->GetDetector() == 57 || | |
291 | trklWord->GetDetector() == 47 || | |
292 | trklWord->GetDetector() == 32) | |
293 | continue; | |
294 | ||
295 | fHistYposRaw->Fill(trklWord->GetY()); | |
296 | fHistdYRaw->Fill(trklWord->GetdY()); | |
297 | fHistZrowRaw->Fill(trklWord->GetZbin()); | |
298 | fHistPidRaw->Fill(trklWord->GetPID()); | |
299 | fHistYdYRaw->Fill(trklWord->GetY(), trklWord->GetdY()); | |
300 | ||
301 | detRaw = trklWord->GetDetector(); | |
302 | robRaw = trklWord->GetROB(); | |
303 | mcmRaw = trklWord->GetMCM(); | |
304 | yRaw = trklWord->GetYbin(); | |
305 | dyRaw = trklWord->GetdY(); | |
306 | trklRaw.Fill(); | |
307 | } | |
308 | } | |
309 | ||
310 | // ----- tracklet comparison raw to simulated ----- | |
311 | trklRaw.BuildIndex("(det+rob*540)*16+mcm", "det"); | |
312 | trklSim.BuildIndex("(det+rob*540)*16+mcm", "det"); | |
313 | trklRaw.AddFriend(&trklSim, (const char*) "sim"); | |
314 | gDirectory->Add(fHistYdiff, kFALSE); | |
315 | Int_t ncomp = trklRaw.Draw("(sim.y-y)*160e-4>>+ydiff", "", "goff"); | |
316 | printf("----- Compared %i tracklets -----\n", ncomp); | |
317 | gDirectory->Remove(fHistYdiff); | |
318 | gDirectory->Add(fHistdYdiff, kFALSE); | |
319 | trklRaw.Draw("(sim.dy-dy)>>+dydiff", "", "goff"); | |
320 | gDirectory->Remove(fHistdYdiff); | |
321 | gDirectory->Add(fHistdYdYraw, kFALSE); | |
322 | trklRaw.Draw("dy:sim.dy>>+dydyraw", "", "goff"); | |
323 | // trklRaw.Scan("det:rob:mcm:y:dy:sim.dy", "sim.dy < 30 && dy > 30"); | |
324 | gDirectory->Remove(fHistdYdYraw); | |
325 | ||
326 | // ----- MC tracks and track references ----- | |
327 | // determine tracklet efficiency | |
328 | if (fMCEvent) { | |
329 | Int_t nTracksMC = fMCEvent->GetNumberOfTracks(); | |
330 | AliInfo(Form("no. of MC tracks: %i", nTracksMC)); | |
331 | for (Int_t iTrack = 0; iTrack < nTracksMC; iTrack++) { | |
332 | // we want primaries | |
333 | if (!fMCEvent->IsPhysicalPrimary(iTrack)) | |
334 | continue; | |
335 | ||
336 | AliMCParticle *mcpart = (AliMCParticle*) fMCEvent->GetTrack(iTrack); | |
337 | ||
338 | // don't look at tracks with too low pt | |
339 | if (TMath::Abs(mcpart->Pt()) < fMinPt) | |
340 | continue; | |
341 | ||
342 | // look for two track references in a chamber | |
343 | Int_t nTrackRefs = mcpart->GetNumberOfTrackReferences(); | |
344 | AliTrackReference *tr[2] = { 0x0 }; | |
345 | Int_t nRef = 0; | |
346 | ||
347 | for (Int_t iTrackRef = 0; iTrackRef < nTrackRefs; iTrackRef++) { | |
348 | AliTrackReference *trackRef = mcpart->GetTrackReference(iTrackRef); | |
349 | if (trackRef->DetectorId() != AliTrackReference::kTRD) | |
350 | continue; | |
351 | if (trackRef->Pt() < fMinPt) | |
352 | continue; | |
353 | Int_t label = trackRef->Label(); | |
354 | if (label < 0) | |
355 | continue; | |
356 | ||
357 | // first track reference, remember it | |
358 | if (nRef == 0) { | |
359 | tr[nRef] = trackRef; | |
360 | nRef++; | |
361 | continue; | |
362 | } | |
363 | else { | |
364 | // next one is too far away, remember it but forget the previous one | |
365 | if (TMath::Abs(trackRef->LocalX() - tr[0]->LocalX()) > 5.) { | |
366 | tr[0] = trackRef; | |
367 | nRef = 1; | |
368 | continue; | |
369 | } | |
370 | // too close to previous track reference | |
371 | // we don't want it | |
372 | else if (TMath::Abs(trackRef->LocalX() - tr[0]->LocalX()) < .5) { | |
373 | continue; | |
374 | } | |
375 | // then it must be ok, so we take it | |
376 | else { | |
377 | tr[1] = trackRef; | |
378 | nRef++; | |
379 | } | |
380 | } | |
381 | ||
382 | // calculation deflection from track references | |
383 | Float_t deflLength = 3. * (tr[1]->LocalY() - tr[0]->LocalY()) / (tr[1]->LocalX() - tr[0]->LocalX()); | |
384 | // if it is too large we reject it | |
385 | if (deflLength < 1.) { | |
386 | fHistCanddY->Fill(deflLength/140e-4); | |
387 | } | |
388 | else { | |
389 | nRef = 0; | |
390 | continue; | |
391 | } | |
392 | ||
393 | // now search for tracklets belonging to this track reference | |
394 | Int_t nTrackletsPerRef = 0; | |
395 | Int_t defl = 0.; | |
396 | for (Int_t iTracklet = 0; iTracklet < fTrackletsSim->GetEntries(); iTracklet++) { | |
397 | trkl = (AliTRDtrackletMCM*) (*fTrackletsSim)[iTracklet]; | |
398 | // they must have the same label | |
399 | if (!trkl->HasLabel(label)) | |
400 | continue; | |
401 | // and be close enough in radial position | |
402 | if (TMath::Abs(trackRef->LocalX() - trkl->GetX()) > 5.) | |
403 | continue; | |
404 | // if they are close in position we accept it | |
405 | if ((TMath::Abs(trackRef->LocalY() - trkl->GetY()) < 5.) && | |
406 | (TMath::Abs(trackRef->Z() - trkl->GetZ()) < 5.)) { | |
407 | defl = trkl->GetdY(); | |
408 | nTrackletsPerRef++; | |
409 | } | |
410 | } | |
411 | fHistTrklPerRef->Fill(nTrackletsPerRef); | |
412 | if (nTrackletsPerRef == 0) { | |
413 | AliInfo(Form("Track ref without assigned tracklet: x=%4.2f, y=%4.2f, z=%4.2f, pt=%4.2f (%i)", | |
414 | trackRef->X(), trackRef->Y(), trackRef->Z(), trackRef->Pt(), trackRef->Label())); | |
415 | } | |
416 | if (nTrackletsPerRef == 1) { | |
417 | fHistdYdYref->Fill(deflLength/140.e-4, defl); | |
418 | fHistFounddY->Fill(deflLength/140.e-4); | |
419 | } | |
420 | nRef = 0; | |
421 | } | |
422 | } | |
423 | } | |
424 | ||
425 | // ----- ESD tracks ----- | |
426 | if (fESD) { | |
427 | for (Int_t iTrack = 0; iTrack < fESD->GetNumberOfTracks(); iTrack++) { | |
428 | AliESDtrack *esdTrack = fESD->GetTrack(iTrack); | |
429 | } | |
430 | } | |
431 | ||
432 | PostData(1, fOutputList); | |
433 | } | |
434 | ||
435 | void AliTRDonlineTrackletQA::LocalInit() | |
436 | { | |
437 | ||
438 | } | |
439 | ||
440 | void AliTRDonlineTrackletQA::Terminate(const Option_t *option) | |
441 | { | |
442 | fOutputList = dynamic_cast<TList*> (GetOutputData(1)); | |
443 | ||
444 | if (!fOutputList) { | |
445 | AliError("No output objects found!"); | |
446 | return; | |
447 | } | |
448 | ||
449 | // fHistYpos = dynamic_cast<TH1F*> (fOutputList->FindObject("ypos")); | |
450 | // if (fHistYpos) { | |
451 | // TCanvas *c = new TCanvas; | |
452 | // c->cd(); | |
453 | // gPad->SetGridx(); | |
454 | // gPad->SetGridy(); | |
455 | // fHistYpos->DrawCopy(); | |
456 | // } | |
457 | } | |
458 | ||
459 | void AliTRDonlineTrackletQA::PlotMC(AliTRDtrackletMCM *trkl) | |
460 | { | |
461 | // compare the given tracklet to the MC information, | |
462 | // i.e. track references | |
463 | ||
464 | Int_t label = trkl->GetLabel(); | |
465 | ||
466 | // check for valid label | |
467 | if (label < 0 ) { | |
468 | AliWarning("MC tracklet has no label"); | |
469 | return; | |
470 | } | |
471 | if (label >= fMCEvent->GetNumberOfTracks()) { | |
472 | AliWarning("MC tracklet has invalid label"); | |
473 | return; | |
474 | } | |
475 | ||
476 | // some cuts on the tracklet | |
477 | // if (!fMCEvent->IsPhysicalPrimary(label)) | |
478 | // return; | |
479 | // if (TMath::Abs(trkl->GetdYdX()) > 0.05) | |
480 | // return; | |
481 | // if (trkl->GetDetector() % 6 != 5) | |
482 | // return; | |
483 | ||
484 | // get MC particle for this tracklet | |
485 | AliMCParticle *track = (AliMCParticle*) fMCEvent->GetTrack(label); | |
486 | ||
487 | // don't look at tracks with too low pt | |
488 | if (TMath::Abs(track->Pt() < fMinPt)) | |
489 | return; | |
490 | ||
491 | // serach track references corresponding to the current tracklet | |
492 | AliTrackReference *tr[2] = { 0x0 }; | |
493 | Int_t nTrackRefs = 0; | |
494 | ||
495 | for (Int_t iTrackRef = 0; iTrackRef < track->GetNumberOfTrackReferences(); iTrackRef++) { | |
496 | AliTrackReference *trackRef = track->GetTrackReference(iTrackRef); | |
497 | if (trackRef->DetectorId() != AliTrackReference::kTRD) | |
498 | continue; | |
499 | if (trackRef->Pt() < fMinPt) | |
500 | continue; | |
501 | ||
502 | if (TMath::Abs(trkl->GetX() - trackRef->LocalX()) > 5.) | |
503 | continue; | |
504 | ||
505 | tr[nTrackRefs++] = trackRef; | |
506 | ||
507 | if (nTrackRefs == 2) | |
508 | break; | |
509 | } | |
510 | ||
511 | // if there were exactly 2 track references | |
512 | // (otherwise something is strange and we want to look at clean cases) | |
513 | // compare tracklet to them | |
514 | if (nTrackRefs == 2) { | |
515 | // sanity check to be in the same sector | |
516 | if ( TMath::Abs((tr[0]->Alpha()*180./TMath::Pi()-10.)/20. - (trkl->GetDetector()/30)) > .1) { | |
517 | AliError("Track reference in different sector"); | |
518 | } | |
519 | // require minimal distance in X and maximum deflection in Y | |
520 | // for the track references | |
521 | else if ((tr[1]->LocalX() - tr[0]->LocalX()) > 0.1 && TMath::Abs(tr[1]->LocalY() - tr[0]->LocalY()) < 1.) { | |
522 | // calculate slope from track references | |
523 | // and check whether it's in the allowed range | |
524 | Float_t slope = 3. * (tr[1]->LocalY() - tr[0]->LocalY()) / (tr[1]->LocalX() - tr[0]->LocalX()); | |
525 | if (TMath::Abs(slope) < 64*140e-4) { | |
526 | AliDebug(1,Form("x1: %f, x0: %f, y1: %f, y0:%f", | |
527 | tr[1]->LocalX(), tr[0]->LocalX(), tr[1]->LocalY(), tr[0]->LocalY() )); | |
528 | // calculate y-position scaled to radial position of the tracklet | |
529 | // and consider the tilting angle of the pads | |
530 | // since the tracklets are affected by it | |
531 | Float_t yMC = (tr[1]->LocalY() + (-0.5+trkl->GetX() - tr[1]->LocalX()) * (tr[1]->LocalY() - tr[0]->LocalY()) / (tr[1]->LocalX() - tr[0]->LocalX())); | |
532 | Float_t yMCtilt = yMC + (TMath::Tan(TMath::Power(-1, (trkl->GetDetector() % 6))*2.*TMath::Pi()/180.) * (tr[1]->Z() - trkl->GetZ())); | |
533 | if (TMath::Abs(trkl->GetY() - yMCtilt) > 10.) { | |
534 | AliError(Form("Deviation too large for tracklet: 0x%08x in det. %i at x = %f, y = %f, z = %f, alpha = %f", | |
535 | trkl->GetTrackletWord(), trkl->GetDetector(), trkl->GetX(), trkl->GetY(), trkl->GetZ(), tr[0]->Alpha())); | |
536 | } | |
537 | fHistYres->Fill(trkl->GetY() - yMCtilt); | |
538 | fHistYresDy->Fill(trkl->GetY() - yMCtilt, trkl->GetdY()); | |
539 | // what about tilt correction here ??? | |
540 | fHistdYres->Fill(3. * trkl->GetdYdX() - | |
541 | 3. * (tr[1]->LocalY() - tr[0]->LocalY()) / (tr[1]->LocalX() - tr[0]->LocalX())); | |
542 | // plot position deviation in pad-coordinates | |
543 | // to study the influence of the position LUT | |
544 | Float_t padWidth = fGeo->GetPadPlane(trkl->GetDetector())->GetWidthIPad(); | |
545 | Float_t yMClocal = yMCtilt/padWidth - floor(yMCtilt/padWidth) - padWidth/2.; | |
546 | Int_t layer = trkl->GetDetector() % 6; | |
547 | fHistYlocal[layer]->Fill(yMClocal, | |
548 | trkl->GetY()/padWidth - floor(trkl->GetY()/padWidth) - padWidth/2. - yMClocal); | |
549 | // and fill everything to the tree | |
550 | fQ0 = trkl->GetQ0(); | |
551 | fQ1 = trkl->GetQ1(); | |
552 | fNHits = trkl->GetNHits(); | |
553 | fYdiff = trkl->GetY() - yMCtilt; | |
554 | fDYdiff = 3. * trkl->GetdYdX() - | |
555 | 3. * (tr[1]->LocalY() - tr[0]->LocalY()) / (tr[1]->LocalX() - tr[0]->LocalX()); | |
556 | fY = trkl->GetY(); | |
557 | fDY = trkl->GetdYdX(); | |
558 | fTreeTracklets->Fill(); | |
559 | // output tracklets with large deviation | |
560 | if (TMath::Abs(fYdiff) > 0.5) { | |
561 | printf("tracklet: y=%4.2f, dy=%4.2f, ydiff=%4.2f, dydiff=%4.2f, q0=%5d, q1=%5d, nhits=%2d, label=%i\n", | |
562 | trkl->GetY(), trkl->GetdYdX(), fYdiff, fDYdiff, fQ0, fQ1, fNHits, label); | |
563 | } | |
564 | } | |
565 | } | |
566 | } | |
567 | } | |
568 | ||
569 | ||
570 | void AliTRDonlineTrackletQA::PlotESD(AliTRDtrackletMCM *trkl) | |
571 | { | |
572 | ||
573 | Float_t xTrkl = trkl->GetX(); | |
574 | Float_t yTrkl = trkl->GetY(); | |
575 | Float_t zTrkl = trkl->GetZ(); | |
576 | ||
577 | Float_t alpha = (trkl->GetDetector() / 30) * 20. + 10.; | |
578 | alpha *= TMath::Pi() / 180.; | |
579 | ||
580 | AliESDEvent *esdEvent = dynamic_cast<AliESDEvent*> (fInputEvent); | |
581 | if (!esdEvent) | |
582 | return; | |
583 | ||
584 | Float_t mag = ((AliMagF*) TGeoGlobalMagField::Instance()->GetField())->SolenoidField(); | |
585 | ||
586 | for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++) { | |
587 | AliESDtrack *track = esdEvent->GetTrack(iTrack); | |
588 | ||
589 | if (!track->GetOuterParam()) | |
590 | continue; | |
591 | ||
592 | AliExternalTrackParam *param = new AliExternalTrackParam(*(track->GetOuterParam())); | |
593 | ||
594 | AliDebug(10, Form("track %i at x = %f", | |
595 | iTrack, param->GetX(), param->GetY())); | |
596 | param->Propagate(alpha, xTrkl, mag); | |
597 | AliDebug(10, Form("after propagating track %i at x = %f, y = %f", | |
598 | iTrack, param->GetX(), param->GetY())); | |
599 | ||
600 | if ((TMath::Abs(xTrkl - param->GetX()) < 10.) && | |
601 | (TMath::Abs(yTrkl - param->GetY()) < 5.) && | |
602 | (TMath::Abs(zTrkl - param->GetZ()) < 10.)) { | |
603 | AliInfo(Form("match of tracklet-track: %i <-> %i", | |
604 | trkl->GetLabel(), track->GetLabel())); | |
605 | AliDebug(5, Form("tracklet position: det: %3i x = %f, y = %f, z = %f, alpha = %f", | |
606 | trkl->GetDetector(), trkl->GetX(), trkl->GetY(), trkl->GetZ(), alpha)); | |
607 | AliDebug(5, Form("after propagating track %i at x = %f, y = %f, z = %f", | |
608 | iTrack, param->GetX(), param->GetY(), param->GetZ())); | |
609 | ||
610 | fHistYresESD->Fill(yTrkl - param->GetY()); | |
611 | } | |
612 | ||
613 | delete param; | |
614 | } | |
615 | ||
616 | } | |
617 | ||
618 | void AliTRDonlineTrackletQA::PlotESD(AliTRDtrackletWord *trkl) | |
619 | { | |
620 | ||
621 | Float_t xTrkl = trkl->GetX(); | |
622 | Float_t yTrkl = trkl->GetY(); | |
623 | Float_t zTrkl = trkl->GetZ(); | |
624 | ||
625 | Float_t alpha = (trkl->GetDetector() / 30) * 20. + 10.; | |
626 | alpha *= TMath::Pi() / 180.; | |
627 | ||
628 | AliESDEvent *esdEvent = dynamic_cast<AliESDEvent*> (fInputEvent); | |
629 | if (!esdEvent) | |
630 | return; | |
631 | ||
632 | Float_t mag = ((AliMagF*) TGeoGlobalMagField::Instance()->GetField())->SolenoidField(); | |
633 | ||
634 | for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++) { | |
635 | AliESDtrack *track = esdEvent->GetTrack(iTrack); | |
636 | ||
637 | if (!track->GetOuterParam()) | |
638 | continue; | |
639 | ||
640 | AliExternalTrackParam *param = new AliExternalTrackParam(*(track->GetOuterParam())); | |
641 | ||
642 | AliDebug(10, Form("track %i at x = %f", | |
643 | iTrack, param->GetX(), param->GetY())); | |
644 | param->Propagate(alpha, xTrkl, mag); | |
645 | AliDebug(10, Form("after propagating track %i at x = %f, y = %f", | |
646 | iTrack, param->GetX(), param->GetY())); | |
647 | ||
648 | if ((TMath::Abs(xTrkl - param->GetX()) < 10.) && | |
649 | (TMath::Abs(yTrkl - param->GetY()) < 5.) && | |
650 | (TMath::Abs(zTrkl - param->GetZ()) < 10.)) { | |
651 | AliDebug(5, Form("tracklet position: det: %3i x = %f, y = %f, z = %f, alpha = %f", | |
652 | trkl->GetDetector(), trkl->GetX(), trkl->GetY(), trkl->GetZ(), alpha)); | |
653 | AliDebug(5, Form("after propagating track %i at x = %f, y = %f, z = %f", | |
654 | iTrack, param->GetX(), param->GetY(), param->GetZ())); | |
655 | ||
656 | fHistYresESD->Fill(yTrkl - param->GetY()); | |
657 | } | |
658 | ||
659 | delete param; | |
660 | } | |
661 | ||
662 | } | |
663 | ||
664 | ||
665 | Int_t AliTRDonlineTrackletQA::GetTrackletsForMC(Int_t label, Int_t idx[]) | |
666 | { | |
667 | return 0; | |
668 | } |