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