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