]>
Commit | Line | Data |
---|---|---|
1 | // $Header$ | |
2 | ||
3 | #include "Track.h" | |
4 | #include "MCHelixLine.hi" | |
5 | #include "PointSet.h" | |
6 | ||
7 | #include <TPolyLine3D.h> | |
8 | #include <TPolyMarker3D.h> | |
9 | #include <TColor.h> | |
10 | ||
11 | // Updates | |
12 | #include <Reve/RGTopFrame.h> | |
13 | #include <TCanvas.h> | |
14 | ||
15 | #include <vector> | |
16 | ||
17 | using namespace Reve; | |
18 | ||
19 | //______________________________________________________________________ | |
20 | // Track | |
21 | // | |
22 | ||
23 | ClassImp(Reve::Track) | |
24 | ||
25 | Track::Track() : | |
26 | Line(), | |
27 | ||
28 | fV(), | |
29 | fP(), | |
30 | fBeta(0), | |
31 | fCharge(0), | |
32 | fLabel(-1), | |
33 | fIndex(-1), | |
34 | fPathMarks(), | |
35 | ||
36 | fRnrStyle(0) | |
37 | {} | |
38 | ||
39 | Track::Track(TParticle* t, Int_t label, TrackRnrStyle* rs): | |
40 | Line(), | |
41 | ||
42 | fV(t->Vx(), t->Vy(), t->Vz()), | |
43 | fP(t->Px(), t->Py(), t->Pz()), | |
44 | fBeta(t->P()/t->Energy()), | |
45 | fCharge(0), | |
46 | fLabel(label), | |
47 | fIndex(-1), | |
48 | fPathMarks(), | |
49 | ||
50 | fRnrStyle(rs) | |
51 | { | |
52 | fLineColor = fRnrStyle->GetColor(); | |
53 | fMainColorPtr = &fLineColor; | |
54 | ||
55 | TParticlePDG* pdgp = t->GetPDG(); | |
56 | if (pdgp) | |
57 | fCharge = (Int_t) TMath::Nint(pdgp->Charge()/3); | |
58 | ||
59 | SetName(t->GetName()); | |
60 | } | |
61 | ||
62 | Track::Track(Reve::MCTrack* t, TrackRnrStyle* rs): | |
63 | Line(), | |
64 | ||
65 | fV(t->Vx(), t->Vy(), t->Vz()), | |
66 | fP(t->Px(), t->Py(), t->Pz()), | |
67 | fBeta(t->P()/t->Energy()), | |
68 | fCharge(0), | |
69 | fLabel(t->label), | |
70 | fIndex(t->index), | |
71 | fPathMarks(), | |
72 | ||
73 | fRnrStyle(rs) | |
74 | { | |
75 | fLineColor = fRnrStyle->GetColor(); | |
76 | fMainColorPtr = &fLineColor; | |
77 | ||
78 | TParticlePDG* pdgp = t->GetPDG(); | |
79 | if(pdgp == 0) { | |
80 | t->ResetPdgCode(); pdgp = t->GetPDG(); | |
81 | } | |
82 | fCharge = (Int_t) TMath::Nint(pdgp->Charge()/3); | |
83 | ||
84 | SetName(t->GetName()); | |
85 | } | |
86 | ||
87 | Track::Track(Reve::RecTrack* t, TrackRnrStyle* rs) : | |
88 | Line(), | |
89 | ||
90 | fV(t->V), | |
91 | fP(t->P), | |
92 | fBeta(t->beta), | |
93 | fCharge(t->sign), | |
94 | fLabel(t->label), | |
95 | fIndex(t->index), | |
96 | fPathMarks(), | |
97 | ||
98 | fRnrStyle(rs) | |
99 | { | |
100 | fLineColor = fRnrStyle->GetColor(); | |
101 | fMainColorPtr = &fLineColor; | |
102 | ||
103 | SetName(t->GetName()); | |
104 | } | |
105 | ||
106 | Track::~Track() | |
107 | { | |
108 | for (vpPathMark_i i=fPathMarks.begin(); i!=fPathMarks.end(); ++i) | |
109 | delete *i; | |
110 | } | |
111 | ||
112 | /* | |
113 | void Track::Reset(Int_t n_points) | |
114 | { | |
115 | delete [] TPolyLine3D::fP; TPolyLine3D::fP = 0; | |
116 | fN = n_points; | |
117 | if(fN) TPolyLine3D::fP = new Float_t [3*fN]; | |
118 | memset(TPolyLine3D::fP, 0, 3*fN*sizeof(Float_t)); | |
119 | fLastPoint = -1; | |
120 | } | |
121 | */ | |
122 | ||
123 | /**************************************************************************/ | |
124 | ||
125 | void Track::MakeTrack() | |
126 | { | |
127 | ||
128 | TrackRnrStyle& RS((fRnrStyle != 0) ? *fRnrStyle : TrackRnrStyle::fgDefStyle); | |
129 | ||
130 | Float_t px = fP.x, py = fP.y, pz = fP.z; | |
131 | ||
132 | MCVertex mc_v0; | |
133 | mc_v0.x = fV.x; | |
134 | mc_v0.y = fV.y; | |
135 | mc_v0.z = fV.z; | |
136 | mc_v0.t = 0; | |
137 | ||
138 | std::vector<MCVertex> track_points; | |
139 | Bool_t decay = kFALSE; | |
140 | ||
141 | if ((TMath::Abs(fV.z) > RS.fMaxZ) || (fV.x*fV.x + fV.y*fV.y > RS.fMaxR*RS.fMaxR)) | |
142 | goto make_polyline; | |
143 | ||
144 | if (fCharge != 0 && TMath::Abs(RS.fMagField) > 1e-5) { | |
145 | ||
146 | // Charged particle in magnetic field | |
147 | ||
148 | Float_t a = RS.fgkB2C * RS.fMagField * fCharge; | |
149 | ||
150 | MCHelix helix(fRnrStyle, &mc_v0, TMath::C()*fBeta, &track_points, a); //m->cm | |
151 | helix.Init(TMath::Sqrt(px*px+py*py), pz); | |
152 | ||
153 | if(!fPathMarks.empty()) | |
154 | { | |
155 | for(std::vector<Reve::PathMark*>::iterator i=fPathMarks.begin(); i!=fPathMarks.end(); ++i) | |
156 | { | |
157 | Reve::PathMark* pm = *i; | |
158 | ||
159 | if (RS.fFitReferences && pm->type == Reve::PathMark::Reference) | |
160 | { | |
161 | if(TMath::Abs(pm->V.z) > RS.fMaxZ | |
162 | || TMath::Sqrt(pm->V.x*pm->V.x + pm->V.y*pm->V.y) > RS.fMaxR ) | |
163 | goto helix_bounds; | |
164 | ||
165 | // printf("%s fit reference \n", fName.Data()); | |
166 | helix.LoopToVertex(px, py, pz, pm->V.x, pm->V.y, pm->V.z); | |
167 | px = pm->P.x; | |
168 | py = pm->P.y; | |
169 | pz = pm->P.z; | |
170 | } | |
171 | else if(RS.fFitDaughters && pm->type == Reve::PathMark::Daughter) | |
172 | { | |
173 | if(TMath::Abs(pm->V.z) > RS.fMaxZ | |
174 | || TMath::Sqrt(pm->V.x*pm->V.x + pm->V.y*pm->V.y) > RS.fMaxR ) | |
175 | goto helix_bounds; | |
176 | ||
177 | // printf("%s fit daughter \n", fName.Data()); | |
178 | helix.LoopToVertex(px, py, pz, pm->V.x, pm->V.y, pm->V.z); | |
179 | px -= pm->P.x; | |
180 | py -= pm->P.y; | |
181 | pz -= pm->P.z; | |
182 | } | |
183 | else if(RS.fFitDecay && pm->type == Reve::PathMark::Decay) | |
184 | { | |
185 | if(TMath::Abs(pm->V.z) > RS.fMaxZ | |
186 | || TMath::Sqrt(pm->V.x*pm->V.x + pm->V.y*pm->V.y) > RS.fMaxR ) | |
187 | goto helix_bounds; | |
188 | helix.LoopToVertex(px, py, pz, pm->V.x, pm->V.y, pm->V.z); | |
189 | decay = true; | |
190 | break; | |
191 | } | |
192 | } | |
193 | } | |
194 | helix_bounds: | |
195 | // go to bounds | |
196 | if(!decay || RS.fFitDecay == kFALSE){ | |
197 | helix.LoopToBounds(px,py,pz); | |
198 | // printf("%s loop to bounds \n",fName.Data() ); | |
199 | } | |
200 | ||
201 | } else { | |
202 | ||
203 | // Neutral particle or no field | |
204 | ||
205 | MCLine line(fRnrStyle, &mc_v0, TMath::C()*fBeta, &track_points); | |
206 | ||
207 | if(!fPathMarks.empty()){ | |
208 | for(std::vector<Reve::PathMark*>::iterator i=fPathMarks.begin(); i!=fPathMarks.end(); ++i) { | |
209 | Reve::PathMark* pm = *i; | |
210 | ||
211 | if(RS.fFitDaughters && pm->type == Reve::PathMark::Daughter){ | |
212 | if(TMath::Abs(pm->V.z) > RS.fMaxZ | |
213 | || TMath::Sqrt(pm->V.x*pm->V.x + pm->V.y*pm->V.y) > RS.fMaxR ) | |
214 | goto line_bounds; | |
215 | line.GotoVertex(pm->V.x, pm->V.y, pm->V.z); | |
216 | fP.x -= pm->P.x; | |
217 | fP.y -= pm->P.y; | |
218 | fP.z -= pm->P.z; | |
219 | } | |
220 | ||
221 | if(RS.fFitDecay && pm->type == Reve::PathMark::Decay){ | |
222 | if(TMath::Abs(pm->V.z) > RS.fMaxZ | |
223 | || TMath::Sqrt(pm->V.x*pm->V.x + pm->V.y*pm->V.y) > RS.fMaxR ) | |
224 | goto line_bounds; | |
225 | line.GotoVertex(pm->V.x, pm->V.y, pm->V.z); | |
226 | decay = true; | |
227 | break; | |
228 | } | |
229 | } | |
230 | } | |
231 | ||
232 | line_bounds: | |
233 | if(!decay || RS.fFitDecay == kFALSE) | |
234 | line.GotoBounds(px,py,pz); | |
235 | ||
236 | } | |
237 | make_polyline: | |
238 | Reset(track_points.size()); | |
239 | for(std::vector<MCVertex>::iterator i=track_points.begin(); i!=track_points.end(); ++i) | |
240 | SetNextPoint(i->x, i->y, i->z); | |
241 | } | |
242 | ||
243 | /**************************************************************************/ | |
244 | ||
245 | void Track::ImportHits() | |
246 | { | |
247 | Reve::LoadMacro("hits_from_label.C"); | |
248 | gROOT->ProcessLine(Form("hits_from_label(%d);", fLabel)); | |
249 | } | |
250 | ||
251 | void Track::ImportClusters() | |
252 | { | |
253 | Reve::LoadMacro("clusters_from_label.C"); | |
254 | gROOT->ProcessLine(Form("clusters_from_label(%d);", fLabel)); | |
255 | } | |
256 | ||
257 | void Track::ImportClustersFromIndex() | |
258 | { | |
259 | static const Exc_t eH("Track::ImportClustersFromIndex "); | |
260 | ||
261 | if (fIndex < 0) | |
262 | throw(eH + "index not set."); | |
263 | ||
264 | Reve::LoadMacro("clusters_from_index.C"); | |
265 | gROOT->ProcessLine(Form("clusters_from_index(%d);", fIndex)); | |
266 | } | |
267 | ||
268 | void Track::ImportDaughters() | |
269 | { | |
270 | static const Exc_t eH("Track::ImportDaughters "); | |
271 | ||
272 | if (fLabel < 0) | |
273 | throw(eH + "label not set."); | |
274 | ||
275 | Reve::LoadMacro("kine_daughter_tracks.C"); | |
276 | gROOT->ProcessLine(Form("kine_daughter_tracks(%d);", fLabel)); | |
277 | } | |
278 | ||
279 | /**************************************************************************/ | |
280 | ||
281 | void Track::PrintKineStack() | |
282 | { | |
283 | Reve::LoadMacro("print_kine_from_label.C"); | |
284 | gROOT->ProcessLine(Form("print_kine_from_label(%d);", fLabel)); | |
285 | } | |
286 | ||
287 | ||
288 | void Track::PrintPathMarks() | |
289 | { | |
290 | static const Exc_t eH("Track::PrintPathMarks "); | |
291 | ||
292 | if (fLabel < 0) | |
293 | throw(eH + "label not set."); | |
294 | ||
295 | printf("Number of path marks %d label %d\n", | |
296 | fPathMarks.size(), fLabel); | |
297 | ||
298 | PathMark* pm; | |
299 | for(vpPathMark_i i=fPathMarks.begin(); i!=fPathMarks.end(); i++) | |
300 | { | |
301 | pm = *i; | |
302 | printf("Reve::PathMark: %-9s p: %8f %8f %8f Vertex: %8e %8e %8e %g \n", | |
303 | pm->type_name(), | |
304 | pm->P.x, pm->P.y, pm->P.z, | |
305 | pm->V.x, pm->V.y, pm->V.z, | |
306 | pm->time); | |
307 | } | |
308 | } | |
309 | ||
310 | /**************************************************************************/ | |
311 | ||
312 | void Track::CtrlClicked(Reve::Track* track) | |
313 | { | |
314 | Emit("CtrlClicked(Reve::Track*)", (Long_t)track); | |
315 | } | |
316 | ||
317 | ||
318 | /**************************************************************************/ | |
319 | /**************************************************************************/ | |
320 | ||
321 | //______________________________________________________________________ | |
322 | // TrackRnrStyle | |
323 | // | |
324 | ||
325 | ClassImp(Reve::TrackRnrStyle) | |
326 | ||
327 | Float_t TrackRnrStyle::fgDefMagField = 5; | |
328 | const Float_t TrackRnrStyle::fgkB2C = 0.299792458e-3; | |
329 | TrackRnrStyle TrackRnrStyle::fgDefStyle; | |
330 | ||
331 | TrackRnrStyle::TrackRnrStyle() : | |
332 | TObject(), | |
333 | ||
334 | fColor(1), | |
335 | fWidth(1), | |
336 | fMagField(fgDefMagField), | |
337 | ||
338 | fMaxR (350), | |
339 | fMaxZ (450), | |
340 | ||
341 | fMaxOrbs (0.5), | |
342 | fMinAng (45), | |
343 | fDelta (0.1), | |
344 | ||
345 | fFitDaughters (kTRUE), | |
346 | fFitReferences (kTRUE), | |
347 | fFitDecay (kTRUE), | |
348 | ||
349 | fRnrDaughters (kTRUE), | |
350 | fRnrReferences (kTRUE), | |
351 | fRnrDecay (kTRUE) | |
352 | {} | |
353 | /**************************************************************************/ | |
354 | /**************************************************************************/ | |
355 | ||
356 | //______________________________________________________________________ | |
357 | // TrackList | |
358 | // | |
359 | ||
360 | ClassImp(Reve::TrackList) | |
361 | ||
362 | void TrackList::Init() | |
363 | { | |
364 | fMarkerStyle = 5; | |
365 | fMarkerColor = 5; | |
366 | // fMarker->SetMarkerSize(0.05); | |
367 | ||
368 | if (fRnrStyle== 0) fRnrStyle = new TrackRnrStyle; | |
369 | SetMainColorPtr(&fRnrStyle->fColor); | |
370 | } | |
371 | ||
372 | TrackList::TrackList(Int_t n_tracks, TrackRnrStyle* rs) : | |
373 | RenderElementListBase(), | |
374 | TPolyMarker3D(n_tracks), | |
375 | ||
376 | fTitle(), | |
377 | ||
378 | fRnrStyle (rs), | |
379 | fRnrTracks (kTRUE), | |
380 | fEditPathMarks (kFALSE) | |
381 | { | |
382 | Init(); | |
383 | } | |
384 | ||
385 | TrackList::TrackList(const Text_t* name, Int_t n_tracks, TrackRnrStyle* rs) : | |
386 | RenderElementListBase(), | |
387 | TPolyMarker3D(n_tracks), | |
388 | ||
389 | fTitle(), | |
390 | ||
391 | fRnrStyle (rs), | |
392 | fRnrTracks (kTRUE) | |
393 | { | |
394 | Init(); | |
395 | SetName(name); | |
396 | } | |
397 | ||
398 | void TrackList::Reset(Int_t n_tracks) | |
399 | { | |
400 | delete [] fP; fP = 0; | |
401 | fN = n_tracks; | |
402 | if(fN) fP = new Float_t [3*fN]; | |
403 | memset(fP, 0, 3*fN*sizeof(Float_t)); | |
404 | fLastPoint = -1; | |
405 | } | |
406 | ||
407 | /**************************************************************************/ | |
408 | ||
409 | void TrackList::Paint(Option_t* option) | |
410 | { | |
411 | if(fRnrElement) { | |
412 | if(fRnrMarkers) { | |
413 | TPolyMarker3D::Paint(option); | |
414 | } | |
415 | if(fRnrTracks) { | |
416 | for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) { | |
417 | if((*i)->GetRnrElement()) | |
418 | (*i)->GetObject()->Paint(option); | |
419 | } | |
420 | } | |
421 | } | |
422 | } | |
423 | ||
424 | /**************************************************************************/ | |
425 | ||
426 | void TrackList::AddElement(RenderElement* el) | |
427 | { | |
428 | static const Exc_t eH("TrackList::AddElement "); | |
429 | if (dynamic_cast<Track*>(el) == 0) | |
430 | throw(eH + "new element not a Track."); | |
431 | RenderElementListBase::AddElement(el); | |
432 | } | |
433 | ||
434 | /**************************************************************************/ | |
435 | ||
436 | void TrackList::MakeTracks() | |
437 | { | |
438 | for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) { | |
439 | ((Track*)(*i))->MakeTrack(); | |
440 | } | |
441 | gReve->Redraw3D(); | |
442 | } | |
443 | ||
444 | ||
445 | void TrackList::MakeMarkers() | |
446 | { | |
447 | Reset(fChildren.size()); | |
448 | for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) { | |
449 | Track& t = *((Track*)(*i)); | |
450 | if(t.GetN() > 0) | |
451 | SetNextPoint(t.fV.x, t.fV.y, t.fV.z); | |
452 | } | |
453 | gReve->Redraw3D(); | |
454 | } | |
455 | ||
456 | /**************************************************************************/ | |
457 | /*************************************************************************/ | |
458 | ||
459 | void TrackList::SetWidth(Width_t w) | |
460 | { | |
461 | Width_t oldw = fRnrStyle->fWidth; | |
462 | fRnrStyle->fWidth = w; | |
463 | for (List_i i=fChildren.begin(); i!=fChildren.end(); ++i) { | |
464 | Track& t = *((Track*)(*i)); | |
465 | if (t.GetLineWidth() == oldw) | |
466 | t.SetLineWidth(w); | |
467 | } | |
468 | } | |
469 | ||
470 | void TrackList::SetMaxR(Float_t x) | |
471 | { | |
472 | fRnrStyle->fMaxR = x; | |
473 | MakeTracks(); | |
474 | MakeMarkers(); | |
475 | } | |
476 | ||
477 | void TrackList::SetMaxZ(Float_t x) | |
478 | { | |
479 | fRnrStyle->fMaxZ = x; | |
480 | MakeTracks(); | |
481 | MakeMarkers(); | |
482 | } | |
483 | ||
484 | void TrackList::SetMaxOrbs(Float_t x) | |
485 | { | |
486 | fRnrStyle->fMaxOrbs = x; | |
487 | MakeTracks(); | |
488 | } | |
489 | ||
490 | void TrackList::SetMinAng(Float_t x) | |
491 | { | |
492 | fRnrStyle->fMinAng = x; | |
493 | MakeTracks(); | |
494 | } | |
495 | ||
496 | void TrackList::SetDelta(Float_t x) | |
497 | { | |
498 | fRnrStyle->fDelta = x; | |
499 | MakeTracks(); | |
500 | } | |
501 | ||
502 | void TrackList::SetFitDaughters(Bool_t x) | |
503 | { | |
504 | fRnrStyle->fFitDaughters = x; | |
505 | MakeTracks(); | |
506 | } | |
507 | ||
508 | void TrackList::SetFitReferences(Bool_t x) | |
509 | { | |
510 | fRnrStyle->fFitReferences = x; | |
511 | MakeTracks(); | |
512 | } | |
513 | ||
514 | void TrackList::SetFitDecay(Bool_t x) | |
515 | { | |
516 | fRnrStyle->fFitDecay = x; | |
517 | MakeTracks(); | |
518 | } | |
519 | ||
520 | void TrackList::SetRnrDecay(Bool_t rnr) | |
521 | { | |
522 | fRnrStyle->fRnrDecay = rnr; | |
523 | MakeTracks(); | |
524 | } | |
525 | ||
526 | void TrackList::SetRnrDaughters(Bool_t rnr) | |
527 | { | |
528 | fRnrStyle->fRnrDaughters = rnr; | |
529 | MakeTracks(); | |
530 | } | |
531 | ||
532 | void TrackList::SetRnrReferences(Bool_t rnr) | |
533 | { | |
534 | fRnrStyle->fRnrReferences = rnr; | |
535 | MakeTracks(); | |
536 | } | |
537 | ||
538 | void TrackList::SetRnrMarkers(Bool_t rnr) | |
539 | { | |
540 | fRnrMarkers = rnr; | |
541 | gReve->Redraw3D(); | |
542 | } | |
543 | ||
544 | void TrackList::SetRnrTracks(Bool_t rnr) | |
545 | { | |
546 | ||
547 | fRnrTracks = rnr; | |
548 | gReve->Redraw3D(); | |
549 | } | |
550 | ||
551 | /**************************************************************************/ | |
552 | /**************************************************************************/ | |
553 | ||
554 | void TrackList::SelectByPt(Float_t min_pt, Float_t max_pt) | |
555 | { | |
556 | Float_t minptsq = min_pt*min_pt; | |
557 | Float_t maxptsq = max_pt*max_pt; | |
558 | Float_t ptsq; | |
559 | ||
560 | for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) { | |
561 | ptsq = ((Track*)(*i))->fP.Perp2(); | |
562 | (*i)->SetRnrElement(ptsq >= minptsq && ptsq <= maxptsq); | |
563 | } | |
564 | } | |
565 | ||
566 | /**************************************************************************/ | |
567 | ||
568 | void TrackList::ImportHits() | |
569 | { | |
570 | for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) { | |
571 | ((Track*)(*i))->ImportHits(); | |
572 | } | |
573 | } | |
574 | ||
575 | void TrackList::ImportClusters() | |
576 | { | |
577 | for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) { | |
578 | ((Track*)(*i))->ImportClusters(); | |
579 | } | |
580 | } | |
581 | ||
582 | /**************************************************************************/ | |
583 | /**************************************************************************/ | |
584 | /**************************************************************************/ | |
585 | ||
586 | #include "RGEditor.h" | |
587 | ||
588 | //______________________________________________________________________ | |
589 | // TrackCounter | |
590 | // | |
591 | ||
592 | ClassImp(TrackCounter) | |
593 | ||
594 | TrackCounter* TrackCounter::fgInstance = 0; | |
595 | ||
596 | TrackCounter::TrackCounter(const Text_t* name, const Text_t* title) : | |
597 | RenderElement(), | |
598 | TNamed(name, title), | |
599 | ||
600 | fBadLineStyle (6), | |
601 | fClickAction (CA_ToggleTrack), | |
602 | fAllTracks (0), | |
603 | fGoodTracks (0), | |
604 | fTrackLists () | |
605 | { | |
606 | if (fgInstance == 0) fgInstance = this; | |
607 | TQObject::Connect("Reve::Track", "CtrlClicked(Reve::Track*)", | |
608 | "Reve::TrackCounter", this, "DoTrackAction(Reve::Track*)"); | |
609 | } | |
610 | ||
611 | TrackCounter::~TrackCounter() | |
612 | { | |
613 | TQObject::Disconnect("Reve::Track", "DoTrackAction(Reve::Track*)"); | |
614 | if (fgInstance == this) fgInstance = 0; | |
615 | } | |
616 | ||
617 | /**************************************************************************/ | |
618 | ||
619 | void TrackCounter::Reset() | |
620 | { | |
621 | printf("TrackCounter::Reset()\n"); | |
622 | fAllTracks = 0; | |
623 | fGoodTracks = 0; | |
624 | TIter next(&fTrackLists); | |
625 | TrackList* tlist; | |
626 | while ((tlist = dynamic_cast<TrackList*>(next()))) | |
627 | tlist->RemoveParent(this); | |
628 | fTrackLists.Clear(); | |
629 | } | |
630 | ||
631 | void TrackCounter::RegisterTracks(TrackList* tlist, Bool_t goodTracks) | |
632 | { | |
633 | // printf("TrackCounter::RegisterTracks '%s', %s\n", | |
634 | // tlist->GetObject()->GetName(), goodTracks ? "good" : "bad"); | |
635 | ||
636 | tlist->AddParent(this); | |
637 | fTrackLists.Add(tlist); | |
638 | ||
639 | List_i i = tlist->BeginChildren(); | |
640 | while (i != tlist->EndChildren()) | |
641 | { | |
642 | Track* t = dynamic_cast<Track*>(*i); | |
643 | if (t != 0) | |
644 | { | |
645 | if (goodTracks) | |
646 | { | |
647 | ++fGoodTracks; | |
648 | } else { | |
649 | t->SetLineStyle(fBadLineStyle); | |
650 | } | |
651 | ++fAllTracks; | |
652 | } | |
653 | ++i; | |
654 | } | |
655 | } | |
656 | ||
657 | void TrackCounter::DoTrackAction(Track* track) | |
658 | { | |
659 | // !!!! No check done if ok. | |
660 | // !!!! Should also override RemoveElementLocal | |
661 | // !!!! But then ... should also sotre local information if track is ok. | |
662 | ||
663 | switch (fClickAction) | |
664 | { | |
665 | ||
666 | case CA_PrintTrackInfo: | |
667 | { | |
668 | printf("Track '%s'\n", track->GetObject()->GetName()); | |
669 | Vector &v = track->fV, &p = track->fP; | |
670 | printf(" Vx=%f, Vy=%f, Vz=%f; Pt=%f, Pz=%f, phi=%f)\n", | |
671 | v.x, v.y, v.z, p.Perp(), p.z, TMath::RadToDeg()*p.Phi()); | |
672 | printf(" <other information should be printed ... full AliESDtrack>\n"); | |
673 | break; | |
674 | } | |
675 | ||
676 | case CA_ToggleTrack: | |
677 | { | |
678 | if (track->GetLineStyle() == 1) | |
679 | { | |
680 | track->SetLineStyle(fBadLineStyle); | |
681 | --fGoodTracks; | |
682 | } else { | |
683 | track->SetLineStyle(1); | |
684 | ++fGoodTracks; | |
685 | } | |
686 | gReve->Redraw3D(); | |
687 | ||
688 | printf("TrackCounter::CountTrack All=%d, Good=%d, Bad=%d\n", | |
689 | fAllTracks, fGoodTracks, fAllTracks-fGoodTracks); | |
690 | ||
691 | if (gReve->GetEditor()->GetModel() == GetObject()) | |
692 | gReve->EditRenderElement(this); | |
693 | ||
694 | break; | |
695 | } | |
696 | ||
697 | } // end switch fClickAction | |
698 | } | |
699 | ||
700 | /**************************************************************************/ | |
701 | ||
702 | void TrackCounter::OutputEventTracks(FILE* out) | |
703 | { | |
704 | if (out == 0) | |
705 | { | |
706 | out = stdout; | |
707 | fprintf(out, "TrackCounter::FinalizeEvent()\n"); | |
708 | } | |
709 | ||
710 | fprintf(out, "Event = %d Ntracks = %d\n", fEventId, fGoodTracks); | |
711 | ||
712 | TIter tlists(&fTrackLists); | |
713 | TrackList* tlist; | |
714 | Int_t cnt = 0; | |
715 | while ((tlist = (TrackList*) tlists()) != 0) | |
716 | { | |
717 | List_i i = tlist->BeginChildren(); | |
718 | while (i != tlist->EndChildren()) | |
719 | { | |
720 | Track* t = dynamic_cast<Track*>(*i); | |
721 | if (t != 0 && t->GetLineStyle() == 1) | |
722 | { | |
723 | ++cnt; | |
724 | fprintf(out, " %2d: chg=%+2d pt=%8.5f eta=%+8.5f\n", | |
725 | cnt, t->fCharge, t->fP.Perp(), t->fP.Eta()); | |
726 | } | |
727 | ++i; | |
728 | } | |
729 | } | |
730 | } |