]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/Reve/Track.cxx
Add getter for fValueIsColor.
[u/mrichter/AliRoot.git] / EVE / Reve / Track.cxx
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( Bool_t recurse)
126 {
127   TrackRnrStyle& RS((fRnrStyle != 0) ? *fRnrStyle : TrackRnrStyle::fgDefStyle);
128
129   Float_t px = fP.x, py = fP.y, pz = fP.z;  
130
131   MCVertex  mc_v0;
132   mc_v0.x = fV.x;
133   mc_v0.y = fV.y; 
134   mc_v0.z = fV.z; 
135   mc_v0.t = 0;
136
137   std::vector<MCVertex> track_points;
138   Bool_t decay = kFALSE;
139
140   if ((TMath::Abs(fV.z) > RS.fMaxZ) || (fV.x*fV.x + fV.y*fV.y > RS.fMaxR*RS.fMaxR)) 
141     goto make_polyline;
142   
143   if (fCharge != 0 && TMath::Abs(RS.fMagField) > 1e-5) {
144
145     // Charged particle in magnetic field
146
147     Float_t a = RS.fgkB2C * RS.fMagField * fCharge;
148    
149     MCHelix helix(fRnrStyle, &mc_v0, TMath::C()*fBeta, &track_points, a); //m->cm
150     helix.Init(TMath::Sqrt(px*px+py*py), pz);
151    
152     if(!fPathMarks.empty())
153     {
154       for(std::vector<Reve::PathMark*>::iterator i=fPathMarks.begin(); i!=fPathMarks.end(); ++i)
155       {
156         Reve::PathMark* pm = *i;
157         
158         if (RS.fFitReferences && pm->type == Reve::PathMark::Reference)
159         {
160           if(TMath::Abs(pm->V.z) > RS.fMaxZ 
161              || TMath::Sqrt(pm->V.x*pm->V.x + pm->V.y*pm->V.y) > RS.fMaxR )
162             goto helix_bounds;
163
164           // printf("%s fit reference  \n", fName.Data()); 
165           helix.LoopToVertex(px, py, pz, pm->V.x, pm->V.y, pm->V.z);
166           px =  pm->P.x;
167           py =  pm->P.y;
168           pz =  pm->P.z;
169         }
170         else if(RS.fFitDaughters &&  pm->type == Reve::PathMark::Daughter)
171         {
172           if(TMath::Abs(pm->V.z) > RS.fMaxZ 
173              || TMath::Sqrt(pm->V.x*pm->V.x + pm->V.y*pm->V.y) > RS.fMaxR )
174             goto helix_bounds;
175
176           // printf("%s fit daughter  \n", fName.Data()); 
177           helix.LoopToVertex(px, py, pz, pm->V.x, pm->V.y, pm->V.z);
178           px -=  pm->P.x;
179           py -=  pm->P.y;
180           pz -=  pm->P.z;
181         }
182         else if(RS.fFitDecay &&  pm->type == Reve::PathMark::Decay)
183         {
184           if(TMath::Abs(pm->V.z) > RS.fMaxZ 
185              || TMath::Sqrt(pm->V.x*pm->V.x + pm->V.y*pm->V.y) > RS.fMaxR )
186             goto helix_bounds;
187           helix.LoopToVertex(px, py, pz, pm->V.x, pm->V.y, pm->V.z);
188           decay = true;
189           break;
190         }
191       }
192     }
193   helix_bounds:
194     // go to bounds
195     if(!decay || RS.fFitDecay == kFALSE){
196       helix.LoopToBounds(px,py,pz);
197       // printf("%s loop to bounds  \n",fName.Data() );
198     }
199
200   } else {
201
202     // Neutral particle or no field
203
204     MCLine line(fRnrStyle, &mc_v0, TMath::C()*fBeta, &track_points);
205    
206     if(!fPathMarks.empty()){
207       for(std::vector<Reve::PathMark*>::iterator i=fPathMarks.begin(); i!=fPathMarks.end(); ++i) {
208         Reve::PathMark* pm = *i;
209
210         if(RS.fFitDaughters &&  pm->type == Reve::PathMark::Daughter){
211           if(TMath::Abs(pm->V.z) > RS.fMaxZ 
212              || TMath::Sqrt(pm->V.x*pm->V.x + pm->V.y*pm->V.y) > RS.fMaxR )
213             goto line_bounds;
214           line.GotoVertex(pm->V.x, pm->V.y, pm->V.z);
215           fP.x -=  pm->P.x;
216           fP.y -=  pm->P.y;
217           fP.z -=  pm->P.z;
218         }
219
220         if(RS.fFitDecay &&  pm->type == Reve::PathMark::Decay){
221           if(TMath::Abs(pm->V.z) > RS.fMaxZ 
222              || TMath::Sqrt(pm->V.x*pm->V.x + pm->V.y*pm->V.y) > RS.fMaxR )
223             goto line_bounds;
224           line.GotoVertex(pm->V.x, pm->V.y, pm->V.z);
225           decay = true;
226           break;
227         }
228       }
229     }
230
231   line_bounds:
232     if(!decay || RS.fFitDecay == kFALSE)
233       line.GotoBounds(px,py,pz);
234
235   }
236 make_polyline:
237   Reset(track_points.size());
238   for(std::vector<MCVertex>::iterator i=track_points.begin(); i!=track_points.end(); ++i)
239     SetNextPoint(i->x, i->y, i->z);
240
241   if(recurse) {
242     Track* t;
243     for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
244       t = dynamic_cast<Track*>(*i);
245       if(t)t->MakeTrack(recurse); 
246     }
247   }
248 }
249
250 /**************************************************************************/
251 namespace {
252 struct cmp_pathmark {
253   bool operator()(PathMark* const & a, PathMark* const & b)
254   { return a->time < b->time; }
255 };
256 }
257
258 void Track::SortPathMarksByTime()
259 {
260  sort(fPathMarks.begin(), fPathMarks.end(), cmp_pathmark());
261 }
262
263
264 /**************************************************************************/
265
266 void Track::ImportHits()
267 {
268   Reve::LoadMacro("hits_from_label.C");
269   gROOT->ProcessLine(Form("hits_from_label(%d, (Reve::RenderElement*)%p);", 
270                           fLabel, this));
271 }
272
273 void Track::ImportClusters()
274 {
275   Reve::LoadMacro("clusters_from_label.C");
276   gROOT->ProcessLine(Form("clusters_from_label(%d, (Reve::RenderElement*)%p);", 
277                           fLabel, this));
278 }
279
280 void Track::ImportClustersFromIndex()
281 {
282   static const Exc_t eH("Track::ImportClustersFromIndex ");
283
284   if (fIndex < 0)
285     throw(eH + "index not set.");
286
287   Reve::LoadMacro("clusters_from_index.C");
288   gROOT->ProcessLine(Form("clusters_from_index(%d, (Reve::RenderElement*)%p);", 
289                           fIndex, this));
290 }
291
292 /**************************************************************************/
293
294 void Track::ImportKine()
295 {
296   static const Exc_t eH("Track::ImportKine ");
297
298   if (fLabel < 0)
299     throw(eH + "label not set.");
300
301   Reve::LoadMacro("kine_tracks.C");
302   gROOT->ProcessLine(Form("kine_track(%d, kFALSE, kTRUE, (Reve::RenderElement*)%p);", 
303                           fLabel, this));
304
305 }
306
307 void Track::ImportKineWithArgs(Bool_t importMother, Bool_t importDaugters)
308 {
309   static const Exc_t eH("Track::ImportKineWithArgs ");
310
311   if (fLabel < 0)
312     throw(eH + "label not set.");
313
314   Reve::LoadMacro("kine_tracks.C");
315   gROOT->ProcessLine(Form("kine_track(%d, %d, %d, (Reve::RenderElement*)%p);", 
316                            fLabel, importMother, importDaugters, this));
317
318 }
319
320 /**************************************************************************/
321
322 void Track::PrintKineStack()
323 {
324   Reve::LoadMacro("print_kine_from_label.C");
325   gROOT->ProcessLine(Form("print_kine_from_label(%d);", fLabel));
326 }
327
328
329 void Track::PrintPathMarks()
330 {
331   static const Exc_t eH("Track::PrintPathMarks ");
332
333   if (fLabel < 0)
334     throw(eH + "label not set.");
335
336   printf("Number of path marks %d label %d\n",
337          fPathMarks.size(), fLabel);
338
339   PathMark* pm;
340   for(vpPathMark_i i=fPathMarks.begin(); i!=fPathMarks.end(); i++) 
341   {
342     pm = *i;
343     printf("Reve::PathMark: %-9s  p: %8f %8f %8f Vertex: %8e %8e %8e %g \n",
344            pm->type_name(),
345            pm->P.x,  pm->P.y, pm->P.z,
346            pm->V.x,  pm->V.y, pm->V.z,
347            pm->time);
348   }
349 }
350
351 /**************************************************************************/
352
353 void Track::CtrlClicked(Reve::Track* track)
354 {
355   Emit("CtrlClicked(Reve::Track*)", (Long_t)track);
356 }
357
358
359 /**************************************************************************/
360 /**************************************************************************/
361
362 //______________________________________________________________________
363 // TrackRnrStyle
364 //
365
366 ClassImp(Reve::TrackRnrStyle)
367
368 Float_t       TrackRnrStyle::fgDefMagField = 5;
369 const Float_t TrackRnrStyle::fgkB2C        = 0.299792458e-3;
370 TrackRnrStyle TrackRnrStyle::fgDefStyle;
371
372 TrackRnrStyle::TrackRnrStyle() :
373   TObject(),
374
375   fColor(1),
376   fWidth(1),
377   fMagField(fgDefMagField),
378
379   fMaxR  (350),
380   fMaxZ  (450),
381
382   fMaxOrbs (0.5),
383   fMinAng  (45),
384   fDelta   (0.1),
385
386   fMinPt   (0.1),
387   fMaxPt   (10),
388
389   fFitDaughters  (kTRUE),
390   fFitReferences (kTRUE),
391   fFitDecay      (kTRUE),
392
393   fRnrDaughters  (kTRUE),
394   fRnrReferences (kTRUE),
395   fRnrDecay      (kTRUE)
396 {}
397 /**************************************************************************/
398 /**************************************************************************/
399
400 //______________________________________________________________________
401 // TrackList
402 //
403
404 ClassImp(Reve::TrackList)
405
406 void TrackList::Init()
407 {
408   fMarkerStyle = 2;
409   fMarkerColor = 4;
410   fMarkerSize  = 0.6;
411
412   if (fRnrStyle== 0) fRnrStyle = new TrackRnrStyle;
413   SetMainColorPtr(&fRnrStyle->fColor);
414 }
415
416 TrackList::TrackList(Int_t n_tracks, TrackRnrStyle* rs) :
417   RenderElement(),
418   TPolyMarker3D(n_tracks),
419
420   fTitle(),
421
422   fRnrStyle      (rs),
423   fRnrTracks     (kTRUE),
424   fEditPathMarks (kFALSE)
425 {
426   Init();
427 }
428
429 TrackList::TrackList(const Text_t* name, Int_t n_tracks, TrackRnrStyle* rs) :
430   RenderElement(),
431   TPolyMarker3D(n_tracks),
432   
433   fTitle(),
434
435   fRnrStyle      (rs),
436   fRnrTracks     (kTRUE),
437   fEditPathMarks (kFALSE)
438 {
439   Init();
440   SetName(name);
441 }
442
443 void TrackList::Reset(Int_t n_tracks)
444 {
445   delete [] fP; fP = 0;
446   fN = n_tracks;
447   if(fN) fP = new Float_t [3*fN];
448   memset(fP, 0, 3*fN*sizeof(Float_t));
449   fLastPoint = -1;
450 }
451
452 /**************************************************************************/
453
454 void TrackList::Paint(Option_t* option)
455 {
456   if(fRnrSelf) {
457     if(fRnrMarkers) {
458       TPolyMarker3D::Paint(option);
459     }
460     if(fRnrTracks && fRnrChildren) {
461       for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
462         if((*i)->GetRnrSelf())
463           (*i)->GetObject()->Paint(option);
464       }
465     }
466   }
467 }
468
469 /**************************************************************************/
470
471 void TrackList::AddElement(RenderElement* el)
472 {
473   static const Exc_t eH("TrackList::AddElement ");
474   if (dynamic_cast<Track*>(el)  == 0)
475     throw(eH + "new element not a Track.");
476   RenderElement::AddElement(el);
477 }
478
479 /**************************************************************************/
480
481 void TrackList::MakeTracks()
482 {
483   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
484     ((Track*)(*i))->MakeTrack();
485   }
486   gReve->Redraw3D();
487 }
488
489
490 void TrackList::MakeMarkers()
491 {
492   Reset(fChildren.size());
493   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
494     Track& t = *((Track*)(*i));
495     if(t.GetN() > 0)
496       SetNextPoint(t.fV.x, t.fV.y, t.fV.z);
497   }
498   gReve->Redraw3D();
499 }
500
501 /**************************************************************************/
502 /*************************************************************************/
503
504 void TrackList::SetWidth(Width_t w)
505 {
506   Width_t oldw = fRnrStyle->fWidth;
507   fRnrStyle->fWidth = w;
508   for (List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
509     Track& t = *((Track*)(*i));
510     if (t.GetLineWidth() == oldw)
511       t.SetLineWidth(w);
512   }
513 }
514
515 void TrackList::SetMaxR(Float_t x)
516 {
517   fRnrStyle->fMaxR = x;
518   MakeTracks();
519   MakeMarkers();
520 }
521
522 void TrackList::SetMaxZ(Float_t x)
523 {
524   fRnrStyle->fMaxZ = x;
525   MakeTracks();
526   MakeMarkers();
527 }
528
529 void TrackList::SetMaxOrbs(Float_t x)
530 {
531   fRnrStyle->fMaxOrbs = x;
532   MakeTracks();
533 }
534
535 void TrackList::SetMinAng(Float_t x)
536 {
537   fRnrStyle->fMinAng = x;
538   MakeTracks();
539 }
540
541 void TrackList::SetDelta(Float_t x)
542 {
543   fRnrStyle->fDelta = x;
544   MakeTracks();
545 }
546
547 void TrackList::SetFitDaughters(Bool_t x)
548 {
549   fRnrStyle->fFitDaughters = x;
550   MakeTracks();
551 }
552
553 void TrackList::SetFitReferences(Bool_t x)
554 {
555   fRnrStyle->fFitReferences = x;
556   MakeTracks();
557 }
558
559 void TrackList::SetFitDecay(Bool_t x)
560 {
561   fRnrStyle->fFitDecay = x;
562   MakeTracks();
563 }
564
565 void TrackList::SetRnrDecay(Bool_t rnr)
566 {
567   fRnrStyle->fRnrDecay = rnr;
568   MakeTracks();
569 }
570
571 void TrackList::SetRnrDaughters(Bool_t rnr)
572 {
573   fRnrStyle->fRnrDaughters = rnr;
574   MakeTracks();
575 }
576
577 void TrackList::SetRnrReferences(Bool_t rnr)
578 {
579   fRnrStyle->fRnrReferences = rnr;
580   MakeTracks();
581 }
582  
583 void TrackList::SetRnrMarkers(Bool_t rnr)
584 {
585   fRnrMarkers = rnr;
586   gReve->Redraw3D();
587 }
588
589 void TrackList::SetRnrTracks(Bool_t rnr)
590 {
591
592   fRnrTracks = rnr;
593   gReve->Redraw3D();
594 }
595
596 /**************************************************************************/
597 /**************************************************************************/
598
599 void TrackList::SelectByPt(Float_t min_pt, Float_t max_pt)
600 {
601   fRnrStyle->fMinPt = min_pt;
602   fRnrStyle->fMaxPt = max_pt;
603
604   Float_t minptsq = min_pt*min_pt;
605   Float_t maxptsq = max_pt*max_pt;
606   Float_t ptsq;
607
608   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
609     ptsq = ((Track*)(*i))->fP.Perp2();
610     (*i)->SetRnrSelf(ptsq >= minptsq && ptsq <= maxptsq);
611   }
612 }
613
614 /**************************************************************************/
615
616 void TrackList::ImportHits()
617 {
618   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
619     ((Track*)(*i))->ImportHits();
620   }
621 }
622
623 void TrackList::ImportClusters()
624 {
625   for(List_i i=fChildren.begin(); i!=fChildren.end(); ++i) {
626     ((Track*)(*i))->ImportClusters();
627   }
628 }
629
630 /**************************************************************************/
631 /**************************************************************************/
632 /**************************************************************************/
633
634 #include "RGEditor.h"
635
636 //______________________________________________________________________
637 // TrackCounter
638 //
639
640 ClassImp(TrackCounter)
641
642 TrackCounter* TrackCounter::fgInstance = 0;
643
644 TrackCounter::TrackCounter(const Text_t* name, const Text_t* title) :
645   RenderElement(),
646   TNamed(name, title),
647
648   fBadLineStyle (6),
649   fClickAction  (CA_ToggleTrack),
650   fAllTracks    (0),
651   fGoodTracks   (0),
652   fTrackLists   ()
653 {
654   if (fgInstance == 0) fgInstance = this;
655   TQObject::Connect("Reve::Track", "CtrlClicked(Reve::Track*)",
656                     "Reve::TrackCounter", this, "DoTrackAction(Reve::Track*)");
657 }
658
659 TrackCounter::~TrackCounter()
660 {
661   TQObject::Disconnect("Reve::Track", "DoTrackAction(Reve::Track*)");
662   if (fgInstance == this) fgInstance = 0;
663 }
664
665 /**************************************************************************/
666
667 void TrackCounter::Reset()
668 {
669   printf("TrackCounter::Reset()\n");
670   fAllTracks  = 0;
671   fGoodTracks = 0;
672   TIter next(&fTrackLists);
673   TrackList* tlist;
674   while ((tlist = dynamic_cast<TrackList*>(next())))
675     tlist->RemoveParent(this);
676   fTrackLists.Clear();
677 }
678
679 void TrackCounter::RegisterTracks(TrackList* tlist, Bool_t goodTracks)
680 {
681   // printf("TrackCounter::RegisterTracks '%s', %s\n",
682   //   tlist->GetObject()->GetName(), goodTracks ? "good" : "bad");
683
684   tlist->AddParent(this);
685   fTrackLists.Add(tlist);
686
687   List_i i = tlist->BeginChildren();
688   while (i != tlist->EndChildren())
689   {
690     Track* t = dynamic_cast<Track*>(*i);
691     if (t != 0)
692     {
693       if (goodTracks)
694       {
695         ++fGoodTracks;
696       } else {
697         t->SetLineStyle(fBadLineStyle);
698       }
699       ++fAllTracks;
700     }
701     ++i;
702   }
703 }
704
705 void TrackCounter::DoTrackAction(Track* track)
706 {
707   // !!!! No check done if ok.
708   // !!!! Should also override RemoveElementLocal
709   // !!!! But then ... should also sotre local information if track is ok.
710
711   switch (fClickAction)
712   {
713
714     case CA_PrintTrackInfo:
715     {
716       printf("Track '%s'\n", track->GetObject()->GetName());
717       Vector &v = track->fV, &p = track->fP;
718       printf("  Vx=%f, Vy=%f, Vz=%f; Pt=%f, Pz=%f, phi=%f)\n",
719              v.x, v.y, v.z, p.Perp(), p.z, TMath::RadToDeg()*p.Phi());
720       printf("  <other information should be printed ... full AliESDtrack>\n");
721       break;
722     }
723
724     case CA_ToggleTrack:
725     {
726       if (track->GetLineStyle() == 1)
727       {
728         track->SetLineStyle(fBadLineStyle);
729         --fGoodTracks;
730       } else {
731         track->SetLineStyle(1);
732         ++fGoodTracks;
733       }
734       gReve->Redraw3D();
735
736       printf("TrackCounter::CountTrack All=%d, Good=%d, Bad=%d\n",
737              fAllTracks, fGoodTracks, fAllTracks-fGoodTracks);
738
739       if (gReve->GetEditor()->GetModel() == GetObject())
740         gReve->EditRenderElement(this);
741
742       break;
743     }
744
745   } // end switch fClickAction
746 }
747
748 /**************************************************************************/
749
750 void TrackCounter::OutputEventTracks(FILE* out)
751 {
752   if (out == 0)
753   {
754     out = stdout;
755     fprintf(out, "TrackCounter::FinalizeEvent()\n");
756   }
757
758   fprintf(out, "Event = %d  Ntracks = %d\n", fEventId, fGoodTracks);
759
760   TIter tlists(&fTrackLists);
761   TrackList* tlist;
762   Int_t cnt = 0;
763   while ((tlist = (TrackList*) tlists()) != 0)
764   {
765     List_i i = tlist->BeginChildren();
766     while (i != tlist->EndChildren())
767     {
768       Track* t = dynamic_cast<Track*>(*i);
769       if (t != 0 && t->GetLineStyle() == 1)
770       {
771         ++cnt;
772         fprintf(out, " %2d: chg=%+2d  pt=%8.5f  eta=%+8.5f\n",
773                cnt, t->fCharge, t->fP.Perp(), t->fP.Eta());
774       }
775       ++i;
776     }
777   }
778 }