]> git.uio.no Git - u/mrichter/AliRoot.git/blame - MUON/AliMUONRecoDisplay.cxx
iModule was defined twice
[u/mrichter/AliRoot.git] / MUON / AliMUONRecoDisplay.cxx
CommitLineData
c7ba256d 1//Authors: Mihaela Gheata, Andrei Gheata 09/10/00
2
3#include <iostream.h>
4#include <AliRun.h>
5#include <TClonesArray.h>
6#include "AliMUONRecoEvent.h"
7#include "AliMUONRecoDisplay.h"
8#include <TROOT.h>
9#include <AliPoints.h>
10#include <TSlider.h>
11#include <TView.h>
12#include <TGeometry.h>
13
14ClassImp(AliMUONRecoDisplay)
15
16//-------------------------------------------------------------------
17AliMUONRecoDisplay::AliMUONRecoDisplay(Int_t nevent)
18 :AliDisplay(750)
19{
20//************ Constructor of the reco. event display**********
21 // get reconstructed event from file
22 fFile = new TFile("tree_reco.root");
23 if (!fFile) {
24 cout << "File tree_reco.root not found\n";
25 gApplication->Terminate(0);
26 }
27 fEvReco = 0;
28 fTree = (TTree *) fFile->Get("TreeRecoEvent");
29 if (!fTree) {
30 cout << "Tree of reconstructed events not found on file. Abort.\n";
31 gApplication->Terminate(0);
32 }
33 if (nevent > fTree->GetEntries()) {
34 cout << "Event number out of range. Aborting.\n";
35 gApplication->Terminate(0);
36 }
37 TBranch *branch = fTree->GetBranch("Event");
38 branch->SetAddress(&fEvReco);
39
40 fEvGen = new AliMUONRecoEvent();
41
42 TFile *galice_file = (TFile*)gROOT->GetListOfFiles()->FindObject("galice.root");
43 galice_file->cd();
44 gAlice->GetEvent(nevent);
45
46
47 fRecoTracks = 0;
48 fGenTracks = 0;
49 fPolyRecoList = 0;
50 fPolyGenList = 0;
51 fHighlited = -1;
52 fMinMomentum = 0;
53 fMaxMomentum = 999;
54 MapEvent(nevent);
55}
56
57//-------------------------------------------------------------------
58AliMUONRecoDisplay::~AliMUONRecoDisplay()
59{
60// Destructor of display object
61 if (fPolyRecoList) {
62 fPolyRecoList->Delete();
63 delete fPolyRecoList;
64 }
65 if (fPolyGenList) {
66 fPolyGenList->Delete();
67 delete fPolyGenList;
68 }
69 delete fEvGen;
70}
71
72//-------------------------------------------------------------------
73void AliMUONRecoDisplay::MapEvent(Int_t nevent)
74{
75// get generated event (nevent) from galice.root and corresponding
76// reconstructed event from tree_reco.root;
77 cout << "mapping event " << nevent << endl;
78 fFile->cd();
79 fTree->GetEntry(nevent);
80 // just testing
81 if (!fEvReco) {
82 cout << "failed to load reco event ! Correct this.\n";
83 gApplication->Terminate(0);
84 }
85 fRecoTracks = fEvReco->TracksPtr();
86 fPolyRecoList = MakePolyLines3D(fRecoTracks);
87 // clear previous event
88 if (fEvGen) fEvGen->Clear();
89 fEvGen->SetNoEvent(nevent);
90 // get list of particles
91 TClonesArray *Particles = gAlice->Particles();
92 // connect MUON module
93 AliDetector *MUON = gAlice->GetDetector("MUON");
94 if (!MUON) {
95 cout << "MUON module not present.\n";
96 gApplication->Terminate(0);
97 }
98 // get the number of generated tracks
99 Int_t ntracks = (Int_t)gAlice->TreeH()->GetEntries();
100 // Fill the fEvGen object
101 AliMUONRecoTrack *gtrack = 0;
102 AliMUONHit *hit = 0;
103 TParticle *particle;
104 Int_t ch;
105 // loop all tracks
106 for (Int_t track=0; track<ntracks; track++) {
107 hit = (AliMUONHit *) MUON->FirstHit(track);
108 if (!hit) continue;
109 particle = (TParticle *) Particles->UncheckedAt(hit->Track());
110 if (IsReconstructible(track) && TMath::Abs(particle->GetPdgCode())==13) {
111 gtrack = fEvGen->AddEmptyTrack();
112 gtrack->SetSign(TMath::Sign((Int_t)1, -particle->GetPdgCode()));
113 // reset hits
114 for (ch=0; ch<10; ch++) gtrack->SetHitPosition(ch,0,0,0);
115 // loop all hits
116 for (AliMUONHit *muonHit=(AliMUONHit*)MUON->FirstHit(track);
117 muonHit;
118 muonHit=(AliMUONHit*)MUON->NextHit()) {
119 ch = muonHit->fChamber - 1;
120 if (ch<0 || ch>9) continue;
121 gtrack->SetHitPosition(ch, muonHit->X(), muonHit->Y(), muonHit->Z());
122 gtrack->SetMomReconstr(particle->Px(), particle->Py(), particle->Pz());
123 gtrack->SetVertexPos(particle->Vz());
124 gtrack->SetChi2r(0.0);
125 }
126 }
127 }
128 fGenTracks = fEvGen->TracksPtr();
129 fPolyGenList = MakePolyLines3D(fGenTracks);
130}
131//-------------------------------------------------------------------
132void AliMUONRecoDisplay::XYPlot()
133{
134// Plot reco. tracks hits in all chambers for current event:
135// - open blue squares : generated muons that were reconstructed accurate
136// - open cyan squares : generated muons that were not reconstructed
137// - filled green circles : reco. tracks (accurate)
138// - filled red squares : fake tracks
139 Double_t kMaxRadius[10];
140 kMaxRadius[0] = kMaxRadius[1] = 91.5;
141 kMaxRadius[2] = kMaxRadius[3] = 122.5;
142 kMaxRadius[4] = kMaxRadius[5] = 158.3;
143 kMaxRadius[6] = kMaxRadius[7] = 260.0;
144 kMaxRadius[8] = kMaxRadius[9] = 260.0;
145
146 TH2F *xygen_found[10];
147 TH2F *xygen_lost[10];
148 TH2F *xyreco_good[10];
149 TH2F *xyreco_fake[10];
150 Double_t x,y,r;
151 Int_t matches[500];
152 Int_t index, ch;
153
154 TPad *pad = (TPad*)gROOT->GetSelectedPad();
155 TCanvas *canvas = new TCanvas("xy", "Reconstruction efficiency");
156 canvas->Clear();
157 canvas->cd();
158 canvas->Divide(3,4);
159 canvas->cd(11);
160 gPad->Delete();
161 canvas->cd(12);
162 gPad->Delete();
163 canvas->cd(1);
164
165 // Define histograms for x-y plots
166 for (ch=0; ch<10; ch++) {
167 xygen_found[ch] = new TH2F("xygen_found","",50,-kMaxRadius[ch],kMaxRadius[ch],50,-kMaxRadius[ch],kMaxRadius[ch]);
168 xygen_lost[ch] = new TH2F("xygen_lost","",50,-kMaxRadius[ch],kMaxRadius[ch],50,-kMaxRadius[ch],kMaxRadius[ch]);
169 xyreco_good[ch] = new TH2F("xyreco_good","",50,-kMaxRadius[ch],kMaxRadius[ch],50,-kMaxRadius[ch],kMaxRadius[ch]);
170 xyreco_fake[ch] = new TH2F("xyreco_fake","",50,-kMaxRadius[ch],kMaxRadius[ch],50,-kMaxRadius[ch],kMaxRadius[ch]);
171 }
172 // find list of matching tracks
173 fPrinted = kTRUE; // no need to print
174 for (index=0; index<fRecoTracks->GetEntriesFast(); index++) {
175 matches[index] = GetBestMatch(index);
176 }
177 // and fill them
178 for (index=0; index<fGenTracks->GetEntriesFast(); index++) { // GEANT
179 Bool_t wasreconst = kFALSE;
180 for (Int_t i=0; i<fRecoTracks->GetEntriesFast(); i++) {
181 if (matches[i] == index) {
182 wasreconst = kTRUE;
183 break;
184 }
185 }
186 AliMUONRecoTrack *current = (AliMUONRecoTrack*)fGenTracks->UncheckedAt(index);
187 for (ch=0; ch<10; ch++) {
188 x = current->GetPosX(ch);
189 y = current->GetPosY(ch);
190 r = TMath::Sqrt(x*x +y*y);
191 if (r >= 10) {
192 if (wasreconst) {
193 xygen_found[ch]->Fill(x,y);
194 } else {xygen_lost[ch]->Fill(x,y);}
195 }
196 }
197 }
198 for (index=0; index<fRecoTracks->GetEntriesFast(); index++) { // reco
199 AliMUONRecoTrack *current = (AliMUONRecoTrack*)fRecoTracks->UncheckedAt(index);
200 for (ch=0; ch<10; ch++) {
201 x = current->GetPosX(ch);
202 y = current->GetPosY(ch);
203 r = TMath::Sqrt(x*x +y*y);
204 if (r >= 10) {
205 if (matches[index] >= 0) {
206 xyreco_good[ch]->Fill(x,y);
207 } else {xyreco_fake[ch]->Fill(x,y);}
208 }
209 }
210 }
211 // finally plot them
212 for (ch=0; ch<10; ch++) {
213 canvas->cd(ch+1);
214 xygen_found[ch]->SetMarkerColor(kBlue);
215 xygen_found[ch]->SetMarkerStyle(4);
216 xygen_found[ch]->SetMarkerSize(0.5);
217 xygen_found[ch]->SetStats(kFALSE);
218 xygen_found[ch]->Draw();
219 xygen_lost[ch]->SetMarkerColor(kCyan);
220 xygen_lost[ch]->SetMarkerStyle(4);
221 xygen_lost[ch]->SetMarkerSize(0.5);
222 xygen_lost[ch]->SetStats(kFALSE);
223 xygen_lost[ch]->Draw("SAME");
224 xyreco_good[ch]->SetMarkerColor(kGreen);
225 xyreco_good[ch]->SetMarkerStyle(20);
226 xyreco_good[ch]->SetMarkerSize(0.4);
227 xyreco_good[ch]->SetStats(kFALSE);
228 xyreco_good[ch]->Draw("SAME");
229 xyreco_fake[ch]->SetMarkerColor(kRed);
230 xyreco_fake[ch]->SetMarkerStyle(20);
231 xyreco_fake[ch]->SetMarkerSize(0.5);
232 xyreco_fake[ch]->SetStats(kFALSE);
233 xyreco_fake[ch]->Draw("SAME");
234 }
235 canvas->SetTitle("y vs. x for simulated and reconstructed tracks");
236 pad->cd();
237}
238
239//-------------------------------------------------------------------
240void AliMUONRecoDisplay::RecoEfficiency(Int_t first, Int_t last)
241{
242// Loop selected reconstructed events, compute efficiency as total number of
243// reconstructed tracks (accurate) over total number of generated
244// reconstructible muons. Make histogram for momentum precision and profile
245// of mean momentum precision versus total momentum (generated)
246 Int_t nevents = (Int_t)fTree->GetEntries();
247 if (last > nevents) last = nevents - 1;
248 if (first < 0) first = 0;
249 nevents = last - first + 1;
250 Int_t track;
251 Float_t efficiency;
252 Float_t generated=0, found=0, fake=0; // number of generated/found/fake tracks
253 Double_t pgen, preco, dP; // dP=preco-pgen
254 AliMUONRecoTrack *rtrack, *gtrack; // generated/reco. current tracks
255 fPrinted = kTRUE; // no need to print
256 Int_t currentEvent = gAlice->GetHeader()->GetEvent();
257
258 TH1F *pReso = new TH1F("Momentum precision", "dP = Prec - Pgen", 100, -5.0, 5.0);
259 pReso->SetXTitle("dP [GeV/c]");
260 pReso->SetYTitle("dN/dP");
261
262 TProfile *dPp = new TProfile("", "dP vs. P", 50, 0, 100, 0, 5);
263 dPp->SetXTitle("P [GeV/c]");
264 dPp->SetYTitle("<dP> [GeV/c]");
265
266 TPad *pad = (TPad*)gROOT->GetSelectedPad();
267 TCanvas *canvas = new TCanvas("c", "Reconstruction efficiency");
268 canvas->Divide(1,2);
269 canvas->cd(1);
270
271 // loop events
272 for (Int_t event=first; event<=last; event++) {
273 // get the reco. & gen. events
274 TFile *galice_file = (TFile*)gROOT->GetListOfFiles()->FindObject("galice.root");
275 galice_file->cd();
276 gAlice->GetEvent(event);
277 MapEvent(event);
278 generated += fGenTracks->GetEntriesFast();
279 // loop reco. tracks
280 for (track=0; track<fRecoTracks->GetEntriesFast(); track++) {
281 rtrack = (AliMUONRecoTrack*)fRecoTracks->UncheckedAt(track);
282 preco = rtrack->P();
283 // find correspondent generated track
284 Int_t ind = GetBestMatch(track);
285 if (ind >=0) {
286 gtrack = (AliMUONRecoTrack*)fGenTracks->UncheckedAt(ind);
287 pgen = gtrack->P();
288 found += 1;
289 dP = preco - pgen;
290 pReso->Fill(dP);
291 dPp->Fill(pgen, TMath::Abs(dP));
292 } else {
293 fake += 1;
294 }
295 }
296 }
297 efficiency = found/generated;
298 cout << "=================================================================\n";
299 cout << "|| Reconstruction efficiency and momentum precision ||\n";
300 cout << "=================================================================\n";
301 cout << "Number of events processed " << nevents << endl;
302 cout << "Total number of reconstructible muons : " << (Int_t)generated << endl;
303 cout << "RECONSTRUCTION EFFICIENCY : " << efficiency*100 << " %" << endl;
304 cout << "Fake track rate : " << 100*fake/generated << " %" << endl;
305 cout << "Momentum precision fit : \n" << endl;
306 pReso->Fit("gaus");
307 cout << "=================================================================\n";
308
309 canvas->cd(1);
310 pReso->Draw();
311 canvas->cd(2);
312 dPp->Draw("HIST");
313 pad->cd();
314 ShowNextEvent(currentEvent-last);
315}
316
317//-------------------------------------------------------------------
318void AliMUONRecoDisplay::ShowNextEvent(Int_t delta)
319{
320// overwritten from AliDisplay in order to get also mapping of next event
321 TFile *galice_file = (TFile*)gROOT->GetListOfFiles()->FindObject("galice.root");
322 galice_file->cd();
323 if (delta) {
324 gAlice->Clear();
325 Int_t currentEvent = gAlice->GetHeader()->GetEvent();
326 Int_t newEvent = currentEvent + delta;
327 Int_t nparticles = gAlice->GetEvent(newEvent);
328 cout << "Event : " << newEvent << " with " << nparticles << " particles\n";
329 if (!gAlice->TreeH()) return;
330 MapEvent(newEvent);
331 fHighlited = -1;
332 }
333 LoadPoints();
334 fPad->cd();
335 Draw();
336 if (gROOT->GetListOfCanvases()->FindObject("xy")) XYPlot();
337}
338//-------------------------------------------------------------------
339Bool_t AliMUONRecoDisplay::IsReconstructible(Int_t track)
340{
341// true if at least three hits in first 2 stations, 3 in last 2 stations
342// and one in station 3
343 AliDetector *MUON = gAlice->GetDetector("MUON");
344 Bool_t chHit[10];
345 Int_t ch;
346 for (ch=0; ch<10; ch++) chHit[ch] = kFALSE;
347 //loop hits
348 for (AliMUONHit *muonHit=(AliMUONHit*)MUON->FirstHit(track);
349 muonHit;
350 muonHit=(AliMUONHit*)MUON->NextHit()) {
351 ch = muonHit->fChamber - 1;
352 if (ch<0 || ch>9) continue;
353 chHit[ch] = kTRUE;
354 }
355 Int_t nhits = 0;
356 for (ch=0; ch<4; ch++) nhits += (chHit[ch])?1:0;
357 if (nhits < 3) return kFALSE;
358 nhits = 0;
359 for (ch=4; ch<6; ch++) nhits+= (chHit[ch])?1:0;
360 if (nhits < 1) return kFALSE;
361
362 for (ch=7; ch<10; ch++) nhits+= (chHit[ch])?1:0;
363 if (nhits < 3) return kFALSE;
364 return kTRUE;
365}
366
367//-------------------------------------------------------------------
368void AliMUONRecoDisplay::DrawView(Float_t theta, Float_t phi, Float_t psi)
369{
370// ovewritten from base class to change the range for MUON
371 gPad->SetCursor(kWatch);
372 gPad->SetFillColor(1);
373 gPad->Clear();
374
375 Int_t iret;
376 TView *view = new TView(1);
377 Float_t range = fRrange*fRangeSlider->GetMaximum();
378 view->SetRange(-range, -range, 0, range, range, 5*range);
379 fZoomX0[0] = -1;
380 fZoomY0[0] = -1;
381 fZoomX1[0] = 1;
382 fZoomY1[0] = 1;
383 // Display Alice geometry
384 gAlice->GetGeometry()->Draw("same");
385 //Loop on all detectors to add their products to the pad
386 DrawHits();
387 // add itself to the list (last)
388 AppendPad();
389
390 view->SetView(phi, theta, psi, iret);
391}
392//-------------------------------------------------------------------
393void AliMUONRecoDisplay::SetDrawHits(Bool_t hits)
394{
395// Turns on/off Geant hits drawing
396 fDrawHits = hits;
397 DrawView(0,-90);
398}
399
400//-------------------------------------------------------------------
401void AliMUONRecoDisplay::CutMomentum(Double_t min, Double_t max)
402{
403// Define momentum cut for displayed reconstructed tracks
404 fMinMomentum = min;
405 fMaxMomentum = max;
406 if (fHighlited >= 0) UnHighlight();
407 DrawView(0,-90);
408}
409
410//-------------------------------------------------------------------
411Int_t AliMUONRecoDisplay::GetBestMatch(Int_t indr, Float_t tolerance)
412{
413// Find the index of best Geant track matching a reconstructed track : track
414// with maximum number of compatible hits (within tolerance*sigma bending and
415// non-bending resolution) and minimum number of fake hits;
416// If no match is found within a given tolerance, the method is called recursively
417// with increasing tolerance, until tolerance = 10;
418 if (indr<0 || indr>=fRecoTracks->GetEntriesFast()) return -1;
419 AliMUONRecoTrack *rtrack = (AliMUONRecoTrack*)fRecoTracks->UncheckedAt(indr);
420 AliMUONRecoTrack *gtrack = 0;
421 Int_t bestMatch = -1;
422 Int_t maxNcompat = 0;
423 Int_t minNfakes = 10;
424 Double_t xrhit,yrhit,radius,xghit,yghit,dX,dY;
425// loop over all Geant tracks
426 for (Int_t indg=0; indg<fGenTracks->GetEntriesFast(); indg++) {
427 gtrack = (AliMUONRecoTrack*)fGenTracks->UncheckedAt(indg);
428 if (!gtrack) continue;
429 Int_t ncompat = 0; // number of compat. hits for this track
430 Int_t nfake = 0; // number of fakes
431 // loop chambers to find compatible hits
432 for (Int_t ch=0; ch<10; ch++) {
433 xrhit = rtrack->GetPosX(ch);
434 yrhit = rtrack->GetPosY(ch);
435 radius = TMath::Sqrt(xrhit*xrhit + yrhit*yrhit);
436 if (radius<10) continue; // skip null hits
437 xghit = gtrack->GetPosX(ch);
438 yghit = gtrack->GetPosY(ch);
439 dX = TMath::Abs(xghit-xrhit);
440 dY = TMath::Abs(yghit-yrhit);
441 if (dX<tolerance*0.144 && dY<tolerance*0.01) { // within 3 sigma resolution
442 ncompat++;
443 continue; // compatible hit
444 } else nfake++; // fake hit
445 }
446 if (ncompat && ncompat>=maxNcompat && nfake<minNfakes) { // this is best matching
447 maxNcompat = ncompat;
448 minNfakes = nfake;
449 bestMatch = indg;
450 }
451 }
452 if (bestMatch<0 && tolerance<=9.) bestMatch = GetBestMatch(indr, tolerance+=1);
453 if (!fPrinted) {
454 rtrack = (AliMUONRecoTrack*)fRecoTracks->UncheckedAt(indr);
455 Int_t sign = rtrack->GetSign();
456 cout << "Reconstructed track : " << indr << "(" << sign << ")" << endl;
457 rtrack->TrackInfo();
458 printf("Best matching Geant track within %i*sgm : %i\n", (Int_t)tolerance, bestMatch);
459 if (bestMatch >=0) {
460 gtrack = (AliMUONRecoTrack*)fGenTracks->UncheckedAt(bestMatch);
461 gtrack->TrackInfo();
462 }
463 cout << "-----------------------------------------------------------------\n";
464 }
465 fPrinted = kTRUE;
466
467 return bestMatch;
468}
469
470//-------------------------------------------------------------------
471void AliMUONRecoDisplay::Highlight(Int_t track)
472{
473// Highlight the specified track
474 if (fHighlited >=0) UnHighlight();
475 if (track<0 || track>fPolyRecoList->GetEntries()) return;
476 TPolyLine3D *line = (TPolyLine3D*)fPolyRecoList->UncheckedAt(track);
477 line->SetLineColor(kYellow);
478 line->SetLineWidth(1);
479 fHighlited = track;
480 DrawView(15,-45,135);
481}
482
483//-------------------------------------------------------------------
484void AliMUONRecoDisplay::UnHighlight()
485{
486// Unhighlight a previous highlighted track
487 if (fHighlited < 0) return; // nothing to do
488 TPolyLine3D *line = (TPolyLine3D*)fPolyRecoList->UncheckedAt(fHighlited);
489 line->SetLineColor(kRed);
490 line->SetLineWidth(1);
491 fHighlited = -1;
492 DrawView(0,-90);
493}
494
495//-------------------------------------------------------------------
496void AliMUONRecoDisplay::DrawHits()
497{
498// Draw hits for all ALICE detectors. Overwrites the DrawHits() method of the
499// base class for reco. track drawing
500
501 Float_t cutmin, cutmax, etamin, etamax, pmom, smin, smax, eta, theta, r;
502 const Float_t kptcutmax = 2;
503 const Float_t ketacutmax = 1.5;
504 Float_t *pxyz;
505 Int_t ntracks,track;
506 TParticle *particle;
507 TObjArray *points;
508 AliPoints *pm;
509
510 //Get cut slider
511 smax = fCutSlider->GetMaximum();
512 smin = fCutSlider->GetMinimum();
513 cutmin = kptcutmax*smin;
514 if (smax < 0.98) cutmax = kptcutmax*smax;
515 else cutmax = 100000;
516
517 //Get eta slider
518 smax = fEtaSlider->GetMaximum();
519 smin = fEtaSlider->GetMinimum();
520 etamin = ketacutmax*(2*smin-1);
521 etamax = ketacutmax*(2*smax-1);
522 if (smin < 0.02) etamin = -1000;
523 if (smax > 0.98) etamax = 1000;
524
525 TIter next(gAlice->Modules());
526 AliModule *module;
527 fHitsCuts = 0;
528 if (fDrawHits) {
529 // draw hits in all modules
530 while((module = (AliModule*)next())) {
531 if (!module->IsActive()) continue;
532 points = module->Points();
533 if (!points) continue;
534 ntracks = points->GetEntriesFast();
535 for (track=0;track<ntracks;track++) {
536 pm = (AliPoints*)points->UncheckedAt(track);
537 if (!pm) continue;
538 particle = pm->GetParticle();
539 if (!particle) continue;
540 pmom = particle->P();
541 if (pmom < cutmin) continue;
542 if (pmom > cutmax) continue;
543 // as a first approximation, take eta of first point
544 pxyz = pm->GetP();
545 r = TMath::Sqrt(pxyz[0]*pxyz[0] + pxyz[1]*pxyz[1]);
546 theta = TMath::ATan2(r,TMath::Abs(pxyz[2]));
547 if(theta) eta = -TMath::Log(TMath::Tan(0.5*theta)); else eta = 1e10;
548 if (pxyz[2] < 0) eta = -eta;
549 if (eta < etamin || eta > etamax) continue;
550 pm->Draw();
551 fHitsCuts += pm->GetN();
552 }
553 }
554 }
555 // draw reconstructed tracks
556 TPolyLine3D *line, *gline;
557 Int_t bestMatch;
558 Double_t px,py,pz,p;
559 AliMUONRecoTrack *rtrack;
560
561 if (fHighlited >= 0) {
562 line = (TPolyLine3D*)fPolyRecoList->UncheckedAt(fHighlited);
563 line->Draw();
564 fPrinted = kFALSE;
565 bestMatch = GetBestMatch(fHighlited);
566 if (bestMatch>=0) {
567 gline = (TPolyLine3D*)fPolyGenList->UncheckedAt(bestMatch);
568 gline->SetLineColor(kRed);
569 gline->SetLineWidth(2);
570 gline->SetLineStyle(2);
571 gline->Draw();
572 }
573 } else {
574 for (track=0; track<fPolyRecoList->GetEntries(); track++) {
575 rtrack = (AliMUONRecoTrack*)fRecoTracks->UncheckedAt(track);
576 px = rtrack->GetMomReconstr(0);
577 py = rtrack->GetMomReconstr(1);
578 pz = rtrack->GetMomReconstr(2);
579 p = rtrack->P();
580 if (p>fMinMomentum && p<fMaxMomentum) {
581 line = (TPolyLine3D*)fPolyRecoList->UncheckedAt(track);
582 line->Draw();
583 }
584 }
585 }
586}
587
588//-------------------------------------------------------------------
589void AliMUONRecoDisplay::ListTracks()
590{
591// List momentum information of all reconstructed traccks within fPmin and fPmax
592// cuts, as well as their best matching Geant tracks
593 cout << "================================================================\n";
594 printf("Reconstructed tracks with momentum in range : %g , %g [GeV/c]\n",
595 fMinMomentum, fMaxMomentum);
596 cout << "----------------------------------------------------------------\n";
597 AliMUONRecoTrack *rtrack;
598 Double_t p;
599 Int_t sign;
600 for (Int_t ind=0; ind<fRecoTracks->GetEntries(); ind++) {
601 rtrack = (AliMUONRecoTrack*)fRecoTracks->UncheckedAt(ind);
602 p = rtrack->P();
603 if (p>fMinMomentum && p<fMaxMomentum) {
604 fPrinted = kFALSE;
605 GetBestMatch(ind);
606 sign = rtrack->GetSign();
607 }
608 }
609 cout << "================================================================\n";
610}
611
612//-------------------------------------------------------------------
613TClonesArray* AliMUONRecoDisplay::MakePolyLines3D(TClonesArray *tracklist)
614{
615// Makes the list of polylines3D corresponding to the list of tracks
616 if (tracklist!=fRecoTracks && tracklist!=fGenTracks) return 0;
617 Bool_t reco = (tracklist==fRecoTracks)?kTRUE:kFALSE;
618 // make sure there is no other list in memory
619 if (reco) {
620 if (fPolyRecoList) {
621 fPolyRecoList->Delete();
622 delete fPolyRecoList;
623 fPolyRecoList = 0;
624 }
625 } else {
626 if (fPolyGenList) {
627 fPolyGenList->Delete();
628 delete fPolyGenList;
629 fPolyGenList = 0;
630 }
631 }
632 if (!tracklist->GetEntries()) return 0;
633
634 AliMUONRecoTrack* track = 0;
635 TClonesArray *polyLines3D = new TClonesArray("TPolyLine3D",1000);
636 TClonesArray &polylist = *polyLines3D;
637 TPolyLine3D *polyline = 0;
638
639 // loop all tracks
640 for (Int_t i=0; i<tracklist->GetEntries(); i++) {
641 track = (AliMUONRecoTrack*)tracklist->UncheckedAt(i);
642 Int_t ch = 0;
643 Double_t x,y,z,r;
644 polyline = new(polylist[i]) TPolyLine3D(2,"");
645 polyline->SetLineColor(kRed);
646 polyline->SetLineWidth(1);
647 polyline->SetNextPoint(0,0,track->GetVertexPos()); // vertex point
648 // loop chambers to fill TPolyLine3D objects
649 for (ch=0; ch<10; ch++) {
650 x = track->GetPosX(ch);
651 y = track->GetPosY(ch);
652 r = TMath::Sqrt(x*x + y*y);
653 if (r < 10) continue;
654 z = track->GetPosZ(ch);
655 polyline->SetNextPoint(x,y,z);
656 }
657 }
658 return polyLines3D;
659}
660
661//-------------------------------------------------------------------
662void AliMUONRecoDisplay::PolyLineInfo(TClonesArray *line3Dlist)
663{
664// Prints information (x, y, z coordinates) for all constructed polylines
665 if (line3Dlist) {
666 TPolyLine3D *polyline = 0;
667 for(Int_t trackIndex=0; trackIndex<line3Dlist->GetEntries(); trackIndex++) {
668 polyline = (TPolyLine3D*)line3Dlist->UncheckedAt(trackIndex);
669 polyline->ls();
670 Float_t *pl = polyline->GetP();
671 for (Int_t i=0; i<polyline->GetN() ;i++) {
672 printf(" x[%d]=%g, y[%d]=%g, z[%d]=%g\n",i,pl[3*i],i,pl[3*i+1],i,pl[3*i+2]);
673 }
674 }
675 }
676}
677
678