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