Logo ROOT   6.18/05
Reference Guide
TEveTrack.cxx
Go to the documentation of this file.
1// @(#)root/eve:$Id$
2// Authors: Matevz Tadel & Alja Mrak-Tadel: 2006, 2007
3
4/*************************************************************************
5 * Copyright (C) 1995-2007, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12#include "TEveTrack.h"
13#include "TEveTrackPropagator.h"
14#include "TEvePointSet.h"
15
16#include "TParticle.h"
17#include "TPolyLine3D.h"
18#include "TMarker.h"
19#include "TPolyMarker3D.h"
20#include "TColor.h"
21#include "TParticlePDG.h"
22
23// Updates
24#include "TEveManager.h"
25#include "TEveBrowser.h"
26#include "TEveTrackProjected.h"
27#include "TEveVSD.h"
28
29#include "Riostream.h"
30
31#include <vector>
32#include <algorithm>
33#include <functional>
34
35/** \class TEveTrack
36\ingroup TEve
37Visual representation of a track.
38
39If member fDpDs is set, the momentum is reduced on all path-marks that do
40not fix the momentum according to the distance travelled from the previous
41pathmark.
42*/
43
45
46////////////////////////////////////////////////////////////////////////////////
47/// Default constructor.
48
50 TEveLine(),
51
52 fV(),
53 fP(),
54 fPEnd(),
55 fBeta(0),
56 fDpDs(0),
57 fPdg(0),
58 fCharge(0),
59 fLabel(kMinInt),
60 fIndex(kMinInt),
61 fStatus(0),
62 fLockPoints(kFALSE),
63 fPathMarks(),
64 fLastPMIdx(0),
65 fPropagator(0)
66{
67}
68
69////////////////////////////////////////////////////////////////////////////////
70/// Constructor from TParticle.
71
73 TEveLine(),
74
75 fV(t->Vx(), t->Vy(), t->Vz()),
76 fP(t->Px(), t->Py(), t->Pz()),
77 fPEnd(),
78 fBeta(t->P()/t->Energy()),
79 fDpDs(0),
80 fPdg(0),
81 fCharge(0),
82 fLabel(label),
83 fIndex(kMinInt),
84 fStatus(t->GetStatusCode()),
85 fLockPoints(kFALSE),
86 fPathMarks(),
87 fLastPMIdx(0),
88 fPropagator(0)
89{
90 SetPropagator(prop);
92
93 TParticlePDG* pdgp = t->GetPDG();
94 if (pdgp) {
95 fPdg = pdgp->PdgCode();
96 fCharge = (Int_t) TMath::Nint(pdgp->Charge()/3);
97 }
98
99 SetName(t->GetName());
100}
101
102////////////////////////////////////////////////////////////////////////////////
103
105 TEveLine(),
106
107 fV(t->Vx(), t->Vy(), t->Vz()),
108 fP(t->Px(), t->Py(), t->Pz()),
109 fPEnd(),
110 fBeta(t->P()/t->Energy()),
111 fDpDs(0),
112 fPdg(0),
113 fCharge(0),
114 fLabel(t->fLabel),
115 fIndex(t->fIndex),
116 fStatus(t->GetStatusCode()),
117 fLockPoints(kFALSE),
118 fPathMarks(),
119 fLastPMIdx(0),
120 fPropagator(0)
121{
122 // Constructor from TEveUtil Monte Carlo track.
123
124 SetPropagator(prop);
126
127 TParticlePDG* pdgp = t->GetPDG();
128 if (pdgp) {
129 fCharge = (Int_t) TMath::Nint(pdgp->Charge()/3);
130 }
131
132 SetName(t->GetName());
133}
134
135////////////////////////////////////////////////////////////////////////////////
136/// Constructor from TEveRecTrack<double> reconstructed track.
137
139 TEveLine(),
140
141 fV(t->fV),
142 fP(t->fP),
143 fPEnd(),
144 fBeta(t->fBeta),
145 fDpDs(0),
146 fPdg(0),
147 fCharge(t->fSign),
148 fLabel(t->fLabel),
149 fIndex(t->fIndex),
150 fStatus(t->fStatus),
151 fLockPoints(kFALSE),
152 fPathMarks(),
153 fLastPMIdx(0),
154 fPropagator(0)
155{
156 SetPropagator(prop);
158
159 SetName(t->GetName());
160}
161
162////////////////////////////////////////////////////////////////////////////////
163/// Constructor from TEveRecTrack<float> reconstructed track.
164/// It is recommended to use constructor with TEveRecTrack<double> since
165/// TEveTrackPropagator operates with double type.
166
168 TEveLine(),
169
170 fV(t->fV),
171 fP(t->fP),
172 fPEnd(),
173 fBeta(t->fBeta),
174 fDpDs(0),
175 fPdg(0),
176 fCharge(t->fSign),
177 fLabel(t->fLabel),
178 fIndex(t->fIndex),
179 fStatus(t->fStatus),
180 fLockPoints(kFALSE),
181 fPathMarks(),
182 fLastPMIdx(0),
183 fPropagator(0)
184{
185 SetPropagator(prop);
187
188 SetName(t->GetName());
189}
190
191////////////////////////////////////////////////////////////////////////////////
192/// Copy constructor. Track parameters are copied but the
193/// extrapolation is not performed so you should still call
194/// MakeTrack() to do that.
195/// If points of 't' are locked, they are cloned.
196
198 TEveLine(),
199 fV(t.fV),
200 fP(t.fP),
201 fPEnd(),
202 fBeta(t.fBeta),
203 fDpDs(t.fDpDs),
204 fPdg(t.fPdg),
205 fCharge(t.fCharge),
206 fLabel(t.fLabel),
207 fIndex(t.fIndex),
208 fStatus(t.fStatus),
209 fLockPoints(t.fLockPoints),
210 fPathMarks(),
211 fLastPMIdx(t.fLastPMIdx),
212 fPropagator(0)
213{
214 if (fLockPoints)
215 ClonePoints(t);
216
217 SetPathMarks(t);
219
220 CopyVizParams(&t);
221}
222
223////////////////////////////////////////////////////////////////////////////////
224/// Destructor.
225
227{
228 SetPropagator(0);
229}
230
231
232////////////////////////////////////////////////////////////////////////////////
233/// Returns list-tree icon for TEveTrack.
234
236{
237 return fgListTreeIcons[4];
238}
239
240////////////////////////////////////////////////////////////////////////////////
241/// Compute the bounding box of the track.
242
244{
245 if (Size() > 0 || ! fPathMarks.empty())
246 {
247 BBoxInit();
248 Int_t n = Size();
250 for (Int_t i = 0; i < n; ++i, p += 3)
251 {
253 }
254 for (vPathMark_ci i = fPathMarks.begin(); i != fPathMarks.end(); ++i)
255 {
256 BBoxCheckPoint(i->fV.fX, i->fV.fY,i->fV.fZ);
257 }
258 }
259 else
260 {
261 BBoxZero();
262 }
263}
264
265////////////////////////////////////////////////////////////////////////////////
266/// Set standard track title based on most data-member values.
267
269{
270 TString idx(fIndex == kMinInt ? "<undef>" : Form("%d", fIndex));
271 TString lbl(fLabel == kMinInt ? "<undef>" : Form("%d", fLabel));
272 SetTitle(Form("Index=%s, Label=%s\nChg=%d, Pdg=%d\n"
273 "pT=%.3f, pZ=%.3f\nV=(%.3f, %.3f, %.3f)",
274 idx.Data(), lbl.Data(), fCharge, fPdg,
275 fP.Perp(), fP.fZ, fV.fX, fV.fY, fV.fZ));
276}
277
278////////////////////////////////////////////////////////////////////////////////
279/// Copy track parameters from t. Track-propagator is set, too.
280/// PathMarks are cleared - you can copy them via SetPathMarks(t).
281/// If track 't' is locked, you should probably clone its points
282/// over - use TEvePointSet::ClonePoints(t);
283
285{
286 fV = t.fV;
287 fP = t.fP;
288 fBeta = t.fBeta;
289 fPdg = t.fPdg;
290 fCharge = t.fCharge;
291 fLabel = t.fLabel;
292 fIndex = t.fIndex;
293
294 fPathMarks.clear();
296}
297
298////////////////////////////////////////////////////////////////////////////////
299/// Copy path-marks from t.
300
302{
303 std::copy(t.RefPathMarks().begin(), t.RefPathMarks().end(),
304 std::back_insert_iterator<vPathMark_t>(fPathMarks));
305}
306
307////////////////////////////////////////////////////////////////////////////////
308/// Set track's render style.
309/// Reference counts of old and new propagator are updated.
310
312{
313 if (fPropagator == prop) return;
315 fPropagator = prop;
317}
318
319////////////////////////////////////////////////////////////////////////////////
320/// Set line and marker attributes from TEveTrackList.
321
323{
324 SetRnrLine(tl->GetRnrLine());
328
333}
334
335////////////////////////////////////////////////////////////////////////////////
336/// Calculate track representation based on track data and current
337/// settings of the propagator.
338/// If recurse is true, descend into children.
339
341{
342 if (!fLockPoints)
343 {
344 Reset(0);
345 fLastPMIdx = 0;
346
348
349 const Double_t maxRsq = rTP.GetMaxR() * rTP.GetMaxR();
350 const Double_t maxZ = rTP.GetMaxZ();
351
352 if ( ! TEveTrackPropagator::IsOutsideBounds(fV, maxRsq, maxZ))
353 {
354 TEveVectorD currP = fP;
355 Bool_t decay = kFALSE;
356 rTP.InitTrack(fV, fCharge);
357 for (vPathMark_i pm = fPathMarks.begin(); pm != fPathMarks.end(); ++pm, ++fLastPMIdx)
358 {
359 Int_t start_point = rTP.GetCurrentPoint();
360
361 if (rTP.GetFitReferences() && pm->fType == TEvePathMarkD::kReference)
362 {
363 if (TEveTrackPropagator::IsOutsideBounds(pm->fV, maxRsq, maxZ))
364 break;
365 if (rTP.GoToVertex(pm->fV, currP))
366 {
367 currP.fX = pm->fP.fX; currP.fY = pm->fP.fY; currP.fZ = pm->fP.fZ;
368 }
369 else
370 {
371 break;
372 }
373 }
374 else if (rTP.GetFitDaughters() && pm->fType == TEvePathMarkD::kDaughter)
375 {
376 if (TEveTrackPropagator::IsOutsideBounds(pm->fV, maxRsq, maxZ))
377 break;
378 if (rTP.GoToVertex(pm->fV, currP))
379 {
380 currP.fX -= pm->fP.fX; currP.fY -= pm->fP.fY; currP.fZ -= pm->fP.fZ;
381 if (fDpDs != 0)
382 {
383 Double_t dp = fDpDs * rTP.GetTrackLength(start_point);
384 Double_t p = currP.Mag();
385 if (p > dp) currP *= 1.0 - dp / p;
386 }
387 }
388 else
389 {
390 break;
391 }
392 }
393 else if (rTP.GetFitDecay() && pm->fType == TEvePathMarkD::kDecay)
394 {
395 if (TEveTrackPropagator::IsOutsideBounds(pm->fV, maxRsq, maxZ))
396 break;
397 rTP.GoToVertex(pm->fV, currP);
398 decay = kTRUE;
399 ++fLastPMIdx;
400 break;
401 }
402 else if (rTP.GetFitCluster2Ds() && pm->fType == TEvePathMarkD::kCluster2D)
403 {
404 TEveVectorD itsect;
405 if (rTP.IntersectPlane(currP, pm->fV, pm->fP, itsect))
406 {
407 TEveVectorD delta = itsect - pm->fV;
408 TEveVectorD vtopass = pm->fV + pm->fE*(pm->fE.Dot(delta));
409 if (TEveTrackPropagator::IsOutsideBounds(vtopass, maxRsq, maxZ))
410 break;
411 if ( ! rTP.GoToVertex(vtopass, currP))
412 break;
413
414 if (fDpDs != 0)
415 {
416 Double_t dp = fDpDs * rTP.GetTrackLength(start_point);
417 Double_t p = currP.Mag();
418 if (p > dp) currP *= 1.0 - dp / p;
419 }
420 }
421 else
422 {
423 Warning("TEveTrack::MakeTrack", "Failed to intersect plane for Cluster2D. Ignoring path-mark.");
424 }
425 }
426 else if (rTP.GetFitLineSegments() && pm->fType == TEvePathMarkD::kLineSegment)
427 {
428 if (TEveTrackPropagator::IsOutsideBounds(pm->fV, maxRsq, maxZ))
429 break;
430
431 if (rTP.GoToLineSegment(pm->fV, pm->fE, currP))
432 {
433 if (fDpDs != 0)
434 {
435 Double_t dp = fDpDs * rTP.GetTrackLength(start_point);
436 Double_t p = currP.Mag();
437 if (p > dp) currP *= 1.0 - dp / p;
438 }
439 }
440 else
441 {
442 break;
443 }
444 }
445 else
446 {
447 if (TEveTrackPropagator::IsOutsideBounds(pm->fV, maxRsq, maxZ))
448 break;
449 }
450 } // loop path-marks
451
452 if (!decay)
453 {
454 // printf("%s loop to bounds \n",fName.Data() );
455 rTP.GoToBounds(currP);
456 }
457 fPEnd = currP;
458 // make_polyline:
459 rTP.FillPointSet(this);
460 rTP.ResetTrack();
461 }
462 }
463
464 if (recurse)
465 {
466 for (List_i i=fChildren.begin(); i!=fChildren.end(); ++i)
467 {
468 TEveTrack* t = dynamic_cast<TEveTrack*>(*i);
469 if (t) t->MakeTrack(recurse);
470 }
471 }
472}
473
474////////////////////////////////////////////////////////////////////////////////
475/// Copy visualization parameters from element el.
476
478{
479 // No local parameters.
480
481 // const TEveTrack* t = dynamic_cast<const TEveTrack*>(el);
482 // if (t)
483 // {}
484
486}
487
488////////////////////////////////////////////////////////////////////////////////
489/// Write visualization parameters.
490
491void TEveTrack::WriteVizParams(std::ostream& out, const TString& var)
492{
493 TEveLine::WriteVizParams(out, var);
494
495 // TString t = " " + var + "->";
496}
497
498////////////////////////////////////////////////////////////////////////////////
499/// Virtual from TEveProjectable, return TEveTrackProjected class.
500
502{
504}
505
506namespace
507{
508 struct Cmp_pathmark_t
509 {
510 bool operator()(TEvePathMarkD const & a, TEvePathMarkD const & b)
511 { return a.fTime < b.fTime; }
512 };
513}
514
515////////////////////////////////////////////////////////////////////////////////
516/// Sort registered pat-marks by time.
517
519{
520 std::sort(fPathMarks.begin(), fPathMarks.end(), Cmp_pathmark_t());
521}
522
523////////////////////////////////////////////////////////////////////////////////
524/// Print registered path-marks.
525
527{
528 static const TEveException eh("TEveTrack::PrintPathMarks ");
529
530 printf("TEveTrack '%s', number of path marks %d, label %d\n",
531 GetName(), (Int_t)fPathMarks.size(), fLabel);
532
533 for (vPathMark_i pm = fPathMarks.begin(); pm != fPathMarks.end(); ++pm)
534 {
535 printf(" %-9s p: %8f %8f %8f Vertex: %8e %8e %8e %g Extra:%8f %8f %8f\n",
536 pm->TypeName(),
537 pm->fP.fX, pm->fP.fY, pm->fP.fZ,
538 pm->fV.fX, pm->fV.fY, pm->fV.fZ,
539 pm->fE.fX, pm->fE.fY, pm->fE.fZ,
540 pm->fTime);
541 }
542}
543
544//------------------------------------------------------------------------------
545
546////////////////////////////////////////////////////////////////////////////////
547/// Emits "SecSelected(TEveTrack*)" signal.
548/// Called from TEveTrackGL on secondary-selection.
549
551{
552 Emit("SecSelected(TEveTrack*)", (Long_t)track);
553}
554
555/** \class TEveTrackList
556\ingroup TEve
557A list of tracks supporting change of common attributes and
558selection based on track parameters.
559*/
560
562
563////////////////////////////////////////////////////////////////////////////////
564/// Constructor. If track-propagator argument is 0, a new default
565/// one is created.
566
569 TAttMarker(1, 20, 1),
570 TAttLine(1,1,1),
571
572 fPropagator(0),
573 fRecurse(kTRUE),
574 fRnrLine(kTRUE),
575 fRnrPoints(kFALSE),
576
577 fMinPt (0), fMaxPt (0), fLimPt (0),
578 fMinP (0), fMaxP (0), fLimP (0)
579{
580
581 fChildClass = TEveTrack::Class(); // override member from base TEveElementList
582
584
585 if (prop == 0) prop = new TEveTrackPropagator;
586 SetPropagator(prop);
587}
588
589////////////////////////////////////////////////////////////////////////////////
590/// Constructor. If track-propagator argument is 0, a new default
591/// one is created.
592
595 TAttMarker(1, 20, 1),
596 TAttLine(1,1,1),
597
598 fPropagator(0),
599 fRecurse(kTRUE),
600 fRnrLine(kTRUE),
601 fRnrPoints(kFALSE),
602
603 fMinPt (0), fMaxPt (0), fLimPt (0),
604 fMinP (0), fMaxP (0), fLimP (0)
605{
606 fChildClass = TEveTrack::Class(); // override member from base TEveElementList
607
609
610 if (prop == 0) prop = new TEveTrackPropagator;
611 SetPropagator(prop);
612}
613
614////////////////////////////////////////////////////////////////////////////////
615/// Destructor.
616
618{
619 SetPropagator(0);
620}
621
622////////////////////////////////////////////////////////////////////////////////
623/// Set default propagator for tracks.
624/// This is not enforced onto the tracks themselves but this is the
625/// propagator that is shown in the TEveTrackListEditor.
626
628{
629 if (fPropagator == prop) return;
631 fPropagator = prop;
633}
634
635////////////////////////////////////////////////////////////////////////////////
636/// Regenerate the visual representations of tracks.
637/// The momentum limits are rescanned during the same traversal.
638
640{
641 fLimPt = fLimP = 0;
642
643 for (List_i i=fChildren.begin(); i!=fChildren.end(); ++i)
644 {
645 TEveTrack* track = dynamic_cast<TEveTrack*>(*i);
646 if (track)
647 {
648 track->MakeTrack(recurse);
649
650 fLimPt = TMath::Max(fLimPt, track->fP.Perp());
651 fLimP = TMath::Max(fLimP, track->fP.Mag());
652 }
653 if (recurse)
654 FindMomentumLimits(*i, recurse);
655 }
656
659
661}
662
663////////////////////////////////////////////////////////////////////////////////
664/// Loop over children and find highest pT and p of contained TEveTracks.
665/// These are stored in members fLimPt and fLimP.
666
668{
669 fLimPt = fLimP = 0;
670
671 if (HasChildren())
672 {
673 for (List_i i=BeginChildren(); i!=EndChildren(); ++i)
674 {
675 TEveTrack* track = dynamic_cast<TEveTrack*>(*i);
676 if (track)
677 {
678 fLimPt = TMath::Max(fLimPt, track->fP.Perp());
679 fLimP = TMath::Max(fLimP, track->fP.Mag());
680 }
681 if (recurse)
682 FindMomentumLimits(*i, recurse);
683 }
684
687 }
688
690}
691
692////////////////////////////////////////////////////////////////////////////////
693/// Loop over track elements of argument el and find highest pT and p.
694/// These are stored in members fLimPt and fLimP.
695
697{
698 for (List_i i=el->BeginChildren(); i!=el->EndChildren(); ++i)
699 {
700 TEveTrack* track = dynamic_cast<TEveTrack*>(*i);
701 if (track)
702 {
703 fLimPt = TMath::Max(fLimPt, track->fP.Perp());
704 fLimP = TMath::Max(fLimP, track->fP.Mag());
705 }
706 if (recurse)
707 FindMomentumLimits(*i, recurse);
708 }
709}
710
711////////////////////////////////////////////////////////////////////////////////
712/// Round the momentum limit up to a nice value.
713
715{
716 using namespace TMath;
717
718 if (x < 1e-3) return 1e-3;
719
720 Double_t fac = Power(10, 1 - Floor(Log10(x)));
721 return Ceil(fac*x) / fac;
722}
723
724////////////////////////////////////////////////////////////////////////////////
725/// Set Min/Max cuts so that they are within detected limits.
726
728{
729 using namespace TMath;
730
732 fMaxPt = fMaxPt == 0 ? fLimPt : Min(fMaxPt, fLimPt);
733 fMinP = Min(fMinP, fLimP);
734 fMaxP = fMaxP == 0 ? fLimP : Min(fMaxP, fLimP);
735}
736
737////////////////////////////////////////////////////////////////////////////////
738/// Set rendering of track as line for the list and the elements.
739
741{
742 for (List_i i=BeginChildren(); i!=EndChildren(); ++i)
743 {
744 TEveTrack* track = (TEveTrack*)(*i);
745 if (track->GetRnrLine() == fRnrLine)
746 track->SetRnrLine(rnr);
747 if (fRecurse)
748 SetRnrLine(rnr, *i);
749 }
750 fRnrLine = rnr;
751}
752
753////////////////////////////////////////////////////////////////////////////////
754/// Set rendering of track as line for children of el.
755
757{
758 TEveTrack* track;
759 for (List_i i=el->BeginChildren(); i!=el->EndChildren(); ++i)
760 {
761 track = dynamic_cast<TEveTrack*>(*i);
762 if (track && (track->GetRnrLine() == fRnrLine))
763 track->SetRnrLine(rnr);
764 if (fRecurse)
765 SetRnrLine(rnr, *i);
766 }
767}
768
769////////////////////////////////////////////////////////////////////////////////
770/// Set rendering of track as points for the list and the elements.
771
773{
774 for (List_i i=BeginChildren(); i!=EndChildren(); ++i)
775 {
776 TEveTrack* track = (TEveTrack*)(*i);
777 if (track->GetRnrPoints() == fRnrPoints)
778 track->SetRnrPoints(rnr);
779 if (fRecurse)
780 SetRnrPoints(rnr, *i);
781 }
782 fRnrPoints = rnr;
783}
784
785////////////////////////////////////////////////////////////////////////////////
786/// Set rendering of track as points for children of el.
787
789{
790 TEveTrack* track;
791 for (List_i i=el->BeginChildren(); i!=el->EndChildren(); ++i)
792 {
793 track = dynamic_cast<TEveTrack*>(*i);
794 if (track)
795 if (track->GetRnrPoints() == fRnrPoints)
796 track->SetRnrPoints(rnr);
797 if (fRecurse)
798 SetRnrPoints(rnr, *i);
799 }
800}
801
802////////////////////////////////////////////////////////////////////////////////
803/// Set main (line) color for the list and the elements.
804
806{
807 for (List_i i=BeginChildren(); i!=EndChildren(); ++i)
808 {
809 TEveTrack* track = (TEveTrack*)(*i);
810 if (track->GetLineColor() == fLineColor)
811 track->SetLineColor(col);
812 if (fRecurse)
813 SetLineColor(col, *i);
814 }
816}
817
818////////////////////////////////////////////////////////////////////////////////
819/// Set line color for children of el.
820
822{
823 TEveTrack* track;
824 for (List_i i=el->BeginChildren(); i!=el->EndChildren(); ++i)
825 {
826 track = dynamic_cast<TEveTrack*>(*i);
827 if (track && track->GetLineColor() == fLineColor)
828 track->SetLineColor(col);
829 if (fRecurse)
830 SetLineColor(col, *i);
831 }
832}
833
834////////////////////////////////////////////////////////////////////////////////
835/// Set line width for the list and the elements.
836
838{
839 for (List_i i=BeginChildren(); i!=EndChildren(); ++i)
840 {
841 TEveTrack* track = (TEveTrack*)(*i);
842 if (track->GetLineWidth() == fLineWidth)
843 track->SetLineWidth(width);
844 if (fRecurse)
845 SetLineWidth(width, *i);
846 }
848}
849
850////////////////////////////////////////////////////////////////////////////////
851/// Set line width for children of el.
852
854{
855 TEveTrack* track;
856 for (List_i i=el->BeginChildren(); i!=el->EndChildren(); ++i)
857 {
858 track = dynamic_cast<TEveTrack*>(*i);
859 if (track && track->GetLineWidth() == fLineWidth)
860 track->SetLineWidth(width);
861 if (fRecurse)
862 SetLineWidth(width, *i);
863 }
864}
865
866////////////////////////////////////////////////////////////////////////////////
867/// Set line style for the list and the elements.
868
870{
871 for (List_i i=BeginChildren(); i!=EndChildren(); ++i)
872 {
873 TEveTrack* track = (TEveTrack*)(*i);
874 if (track->GetLineStyle() == fLineStyle)
875 track->SetLineStyle(style);
876 if (fRecurse)
877 SetLineStyle(style, *i);
878 }
880}
881
882////////////////////////////////////////////////////////////////////////////////
883/// Set line style for children of el.
884
886{
887 TEveTrack* track;
888 for (List_i i=el->BeginChildren(); i!=el->EndChildren(); ++i)
889 {
890 track = dynamic_cast<TEveTrack*>(*i);
891 if (track && track->GetLineStyle() == fLineStyle)
892 track->SetLineStyle(style);
893 if (fRecurse)
894 SetLineStyle(style, *i);
895 }
896}
897
898////////////////////////////////////////////////////////////////////////////////
899/// Set marker style for the list and the elements.
900
902{
903 for (List_i i=BeginChildren(); i!=EndChildren(); ++i)
904 {
905 TEveTrack* track = (TEveTrack*)(*i);
906 if (track->GetMarkerStyle() == fMarkerStyle)
907 track->SetMarkerStyle(style);
908 if (fRecurse)
910 }
912}
913
914////////////////////////////////////////////////////////////////////////////////
915/// Set marker style for children of el.
916
918{
919 TEveTrack* track;
920 for (List_i i=el->BeginChildren(); i!=el->EndChildren(); ++i)
921 {
922 track = dynamic_cast<TEveTrack*>(*i);
923 if (track && track->GetMarkerStyle() == fMarkerStyle)
924 track->SetMarkerStyle(style);
925 if (fRecurse)
927 }
928}
929
930////////////////////////////////////////////////////////////////////////////////
931/// Set marker color for the list and the elements.
932
934{
935 for (List_i i=BeginChildren(); i!=EndChildren(); ++i)
936 {
937 TEveTrack* track = (TEveTrack*)(*i);
938 if (track->GetMarkerColor() == fMarkerColor)
939 track->SetMarkerColor(col);
940 if (fRecurse)
941 SetMarkerColor(col, *i);
942 }
943 fMarkerColor = col;
944}
945
946////////////////////////////////////////////////////////////////////////////////
947/// Set marker color for children of el.
948
950{
951 TEveTrack* track;
952 for (List_i i=el->BeginChildren(); i!=el->EndChildren(); ++i)
953 {
954 track = dynamic_cast<TEveTrack*>(*i);
955 if (track && track->GetMarkerColor() == fMarkerColor)
956 track->SetMarkerColor(col);
957 if (fRecurse)
958 SetMarkerColor(col, *i);
959 }
960}
961
962////////////////////////////////////////////////////////////////////////////////
963/// Set marker size for the list and the elements.
964
966{
967 for (List_i i=BeginChildren(); i!=EndChildren(); ++i)
968 {
969 TEveTrack* track = (TEveTrack*)(*i);
970 if (track->GetMarkerSize() == fMarkerSize)
971 track->SetMarkerSize(size);
972 if (fRecurse)
973 SetMarkerSize(size, *i);
974 }
975 fMarkerSize = size;
976}
977
978////////////////////////////////////////////////////////////////////////////////
979/// Set marker size for children of el.
980
982{
983 TEveTrack* track;
984 for (List_i i=el->BeginChildren(); i!=el->EndChildren(); ++i)
985 {
986 track = dynamic_cast<TEveTrack*>(*i);
987 if (track && track->GetMarkerSize() == fMarkerSize)
988 track->SetMarkerSize(size);
989 if (fRecurse)
990 SetMarkerSize(size, *i);
991 }
992}
993
994////////////////////////////////////////////////////////////////////////////////
995/// Select visibility of tracks by transverse momentum.
996/// If data-member fRecurse is set, the selection is applied
997/// recursively to all children.
998
1000{
1001 fMinPt = min_pt;
1002 fMaxPt = max_pt;
1003
1004 const Double_t minptsq = min_pt*min_pt;
1005 const Double_t maxptsq = max_pt*max_pt;
1006
1007 for (List_i i=BeginChildren(); i!=EndChildren(); ++i)
1008 {
1009 const Double_t ptsq = ((TEveTrack*)(*i))->fP.Perp2();
1010 Bool_t on = ptsq >= minptsq && ptsq <= maxptsq;
1011 (*i)->SetRnrState(on);
1012 if (on && fRecurse)
1013 SelectByPt(min_pt, max_pt, *i);
1014 }
1015}
1016
1017////////////////////////////////////////////////////////////////////////////////
1018/// Select visibility of el's children tracks by transverse momentum.
1019
1021{
1022 const Double_t minptsq = min_pt*min_pt;
1023 const Double_t maxptsq = max_pt*max_pt;
1024
1025 for (List_i i=el->BeginChildren(); i!=el->EndChildren(); ++i)
1026 {
1027 TEveTrack* track = dynamic_cast<TEveTrack*>(*i);
1028 if (track)
1029 {
1030 const Double_t ptsq = track->fP.Perp2();
1031 Bool_t on = ptsq >= minptsq && ptsq <= maxptsq;
1032 track->SetRnrState(on);
1033 if (on && fRecurse)
1034 SelectByPt(min_pt, max_pt, *i);
1035 }
1036 }
1037}
1038
1039////////////////////////////////////////////////////////////////////////////////
1040/// Select visibility of tracks by momentum.
1041/// If data-member fRecurse is set, the selection is applied
1042/// recursively to all children.
1043
1045{
1046 fMinP = min_p;
1047 fMaxP = max_p;
1048
1049 const Double_t minpsq = min_p*min_p;
1050 const Double_t maxpsq = max_p*max_p;
1051
1052 for (List_i i=BeginChildren(); i!=EndChildren(); ++i)
1053 {
1054 const Double_t psq = ((TEveTrack*)(*i))->fP.Mag2();
1055 Bool_t on = psq >= minpsq && psq <= maxpsq;
1056 (*i)->SetRnrState(psq >= minpsq && psq <= maxpsq);
1057 if (on && fRecurse)
1058 SelectByP(min_p, max_p, *i);
1059 }
1060}
1061
1062////////////////////////////////////////////////////////////////////////////////
1063/// Select visibility of el's children tracks by momentum.
1064
1066{
1067 const Double_t minpsq = min_p*min_p;
1068 const Double_t maxpsq = max_p*max_p;
1069
1070 for (List_i i=el->BeginChildren(); i!=el->EndChildren(); ++i)
1071 {
1072 TEveTrack* track = dynamic_cast<TEveTrack*>(*i);
1073 if (track)
1074 {
1075 const Double_t psq = ((TEveTrack*)(*i))->fP.Mag2();
1076 Bool_t on = psq >= minpsq && psq <= maxpsq;
1077 track->SetRnrState(on);
1078 if (on && fRecurse)
1079 SelectByP(min_p, max_p, *i);
1080 }
1081 }
1082}
1083
1084////////////////////////////////////////////////////////////////////////////////
1085/// Find track by label, select it and display it in the editor.
1086
1088{
1089 for (List_i i=fChildren.begin(); i!=fChildren.end(); ++i)
1090 {
1091 if (((TEveTrack*)(*i))->GetLabel() == label)
1092 {
1094 TGListTreeItem *mlti = lt->GetSelected();
1095 if (mlti->GetUserData() != this)
1096 mlti = FindListTreeItem(lt);
1097 TGListTreeItem *tlti = (*i)->FindListTreeItem(lt, mlti);
1098 lt->HighlightItem(tlti);
1099 lt->SetSelected(tlti);
1100 gEve->EditElement(*i);
1101 return (TEveTrack*) *i;
1102 }
1103 }
1104 return 0;
1105}
1106
1107////////////////////////////////////////////////////////////////////////////////
1108/// Find track by index, select it and display it in the editor.
1109
1111{
1112 for (List_i i=fChildren.begin(); i!=fChildren.end(); ++i)
1113 {
1114 if (((TEveTrack*)(*i))->GetIndex() == index)
1115 {
1117 TGListTreeItem *mlti = lt->GetSelected();
1118 if (mlti->GetUserData() != this)
1119 mlti = FindListTreeItem(lt);
1120 TGListTreeItem *tlti = (*i)->FindListTreeItem(lt, mlti);
1121 lt->HighlightItem(tlti);
1122 lt->SetSelected(tlti);
1123 gEve->EditElement(*i);
1124 return (TEveTrack*) *i;
1125 }
1126 }
1127 return 0;
1128}
1129
1130////////////////////////////////////////////////////////////////////////////////
1131/// Copy visualization parameters from element el.
1132
1134{
1135 const TEveTrackList* m = dynamic_cast<const TEveTrackList*>(el);
1136 if (m)
1137 {
1140 fRecurse = m->fRecurse;
1141 fRnrLine = m->fRnrLine;
1142 fRnrPoints = m->fRnrPoints;
1143 fMinPt = m->fMinPt;
1144 fMaxPt = m->fMaxPt;
1145 fLimPt = m->fLimPt;
1146 fMinP = m->fMinP;
1147 fMaxP = m->fMaxP;
1148 fLimP = m->fLimP;
1149 }
1150
1152}
1153
1154////////////////////////////////////////////////////////////////////////////////
1155/// Write visualization parameters.
1156
1157void TEveTrackList::WriteVizParams(std::ostream& out, const TString& var)
1158{
1160
1161 TString t = " " + var + "->";
1164 out << t << "SetRecurse(" << ToString(fRecurse) << ");\n";
1165 out << t << "SetRnrLine(" << ToString(fRnrLine) << ");\n";
1166 out << t << "SetRnrPoints(" << ToString(fRnrPoints) << ");\n";
1167 // These setters are not available -- need proper AND/OR mode.
1168 // out << t << "SetMinPt(" << fMinPt << ");\n";
1169 // out << t << "SetMaxPt(" << fMaxPt << ");\n";
1170 // out << t << "SetLimPt(" << fLimPt << ");\n";
1171 // out << t << "SetMinP(" << fMinP << ");\n";
1172 // out << t << "SetMaxP(" << fMaxP << ");\n";
1173 // out << t << "SetLimP(" << fLimP << ");\n";
1174}
1175
1176////////////////////////////////////////////////////////////////////////////////
1177/// Virtual from TEveProjectable, returns TEveTrackListProjected class.
1178
1180{
1182}
void Class()
Definition: Class.C:29
#define b(i)
Definition: RSha256.hxx:100
#define e(i)
Definition: RSha256.hxx:103
const Int_t kMinInt
Definition: RtypesCore.h:100
int Int_t
Definition: RtypesCore.h:41
float Size_t
Definition: RtypesCore.h:83
const Bool_t kFALSE
Definition: RtypesCore.h:88
long Long_t
Definition: RtypesCore.h:50
short Width_t
Definition: RtypesCore.h:78
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
short Color_t
Definition: RtypesCore.h:79
short Style_t
Definition: RtypesCore.h:76
float Float_t
Definition: RtypesCore.h:53
const Bool_t kTRUE
Definition: RtypesCore.h:87
#define ClassImp(name)
Definition: Rtypes.h:365
include TDocParser_001 C image html pict1_TDocParser_001 png width
Definition: TDocParser.cxx:121
R__EXTERN TEveManager * gEve
Definition: TEveManager.h:243
char name[80]
Definition: TGX11.cxx:109
TRObject operator()(const T1 &t1) const
Binding & operator=(OUT(*fun)(void))
char * Form(const char *fmt,...)
void BBoxCheckPoint(Float_t x, Float_t y, Float_t z)
Definition: TAttBBox.h:58
void BBoxZero(Float_t epsilon=0, Float_t x=0, Float_t y=0, Float_t z=0)
Create cube of volume (2*epsilon)^3 at (x,y,z).
Definition: TAttBBox.cxx:42
void BBoxInit(Float_t infinity=1e6)
Dynamic Float_t[6] X(min,max), Y(min,max), Z(min,max)
Definition: TAttBBox.cxx:29
Line Attributes class.
Definition: TAttLine.h:18
virtual Color_t GetLineColor() const
Return the line color.
Definition: TAttLine.h:33
virtual Width_t GetLineWidth() const
Return the line width.
Definition: TAttLine.h:35
Width_t fLineWidth
Line width.
Definition: TAttLine.h:23
virtual Style_t GetLineStyle() const
Return the line style.
Definition: TAttLine.h:34
Style_t fLineStyle
Line style.
Definition: TAttLine.h:22
Color_t fLineColor
Line color.
Definition: TAttLine.h:21
virtual void SaveLineAttributes(std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1, Int_t widdef=1)
Save line attributes as C++ statement(s) on output stream out.
Definition: TAttLine.cxx:262
Marker Attributes class.
Definition: TAttMarker.h:19
virtual void SaveMarkerAttributes(std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1, Int_t sizdef=1)
Save line attributes as C++ statement(s) on output stream out.
Definition: TAttMarker.cxx:245
virtual Style_t GetMarkerStyle() const
Return the marker style.
Definition: TAttMarker.h:32
virtual Color_t GetMarkerColor() const
Return the marker color.
Definition: TAttMarker.h:31
Color_t fMarkerColor
Marker color.
Definition: TAttMarker.h:22
virtual Size_t GetMarkerSize() const
Return the marker size.
Definition: TAttMarker.h:33
Size_t fMarkerSize
Marker size.
Definition: TAttMarker.h:24
Style_t fMarkerStyle
Marker style.
Definition: TAttMarker.h:23
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition: TClass.h:75
A list of TEveElements.
Definition: TEveElement.h:431
TClass * fChildClass
Definition: TEveElement.h:437
Base class for TEveUtil visualization elements, providing hierarchy management, rendering control and...
Definition: TEveElement.h:34
List_t fChildren
Definition: TEveElement.h:79
List_i EndChildren()
Definition: TEveElement.h:165
Bool_t HasChildren() const
Definition: TEveElement.h:169
virtual TGListTreeItem * FindListTreeItem(TGListTree *ltree)
Find any list-tree-item of this element in list-tree 'ltree'.
virtual void SetMainColor(Color_t color)
Set main color of the element.
Color_t * fMainColorPtr
Definition: TEveElement.h:97
virtual Bool_t SetRnrState(Bool_t rnr)
Set render state of this element and of its children to the same value.
virtual void CopyVizParams(const TEveElement *el)
Copy visualization parameters from element el.
static const char * ToString(Bool_t b)
Convert Bool_t to string - kTRUE or kFALSE.
static const TGPicture * fgListTreeIcons[9]
Definition: TEveElement.h:63
virtual void WriteVizParams(std::ostream &out, const TString &var)
Write-out visual parameters for this object.
List_i BeginChildren()
Definition: TEveElement.h:164
List_t::iterator List_i
Definition: TEveElement.h:70
Exception class thrown by TEve classes and macros.
Definition: TEveUtil.h:103
TGListTree * GetListTree() const
Definition: TEveBrowser.h:114
An arbitrary polyline with fixed line and marker attributes.
Definition: TEveLine.h:26
virtual void SetLineStyle(Style_t lstyle)
Set line-style of the line.
Definition: TEveLine.cxx:86
virtual void SetLineWidth(Width_t lwidth)
Set line-style of the line.
Definition: TEveLine.cxx:106
virtual void SetLineColor(Color_t col)
Set the line color.
Definition: TEveLine.h:48
Bool_t GetRnrPoints() const
Definition: TEveLine.h:53
void SetRnrLine(Bool_t r)
Set rendering of line. Propagate to projected lines.
Definition: TEveLine.cxx:125
Bool_t GetRnrLine() const
Definition: TEveLine.h:52
virtual void CopyVizParams(const TEveElement *el)
Copy visualization parameters from element el.
Definition: TEveLine.cxx:260
virtual void WriteVizParams(std::ostream &out, const TString &var)
Write visualization parameters.
Definition: TEveLine.cxx:277
void SetRnrPoints(Bool_t r)
Set rendering of points. Propagate to projected lines.
Definition: TEveLine.cxx:144
virtual void SetMarkerColor(Color_t col)
Set marker color. Propagate to projected lines.
Definition: TEveLine.cxx:66
TEveGListTreeEditorFrame * GetLTEFrame() const
Definition: TEveManager.h:138
void EditElement(TEveElement *element)
Show element in default editor.
Special-point on track:
Definition: TEvePathMark.h:23
virtual void SetMarkerStyle(Style_t mstyle=1)
Set marker style, propagate to projecteds.
virtual void ClonePoints(const TEvePointSet &e)
Clone points and all point-related information from point-set 'e'.
virtual void SetTitle(const char *t)
Definition: TEvePointSet.h:69
void Reset(Int_t n_points=0, Int_t n_int_ids=0)
Drop all data and set-up the data structures to recive new data.
virtual void SetMarkerSize(Size_t msize=1)
Set marker size, propagate to projecteds.
Base-class for non-linear projections.
virtual void IncRefCount(TEveElement *re)
Increase reference count and add re to the list of back-references.
Definition: TEveUtil.cxx:581
virtual void DecRefCount(TEveElement *re)
Decrease reference count and remove re from the list of back-references.
Definition: TEveUtil.cxx:590
A list of tracks supporting change of common attributes and selection based on track parameters.
Definition: TEveTrack.h:140
virtual void WriteVizParams(std::ostream &out, const TString &var)
Write visualization parameters.
Definition: TEveTrack.cxx:1157
virtual TClass * ProjectedClass(const TEveProjection *p) const
Virtual from TEveProjectable, returns TEveTrackListProjected class.
Definition: TEveTrack.cxx:1179
Double_t fMinPt
Definition: TEveTrack.h:155
virtual void SetLineWidth(Width_t w)
Set line width for the list and the elements.
Definition: TEveTrack.cxx:837
void SetPropagator(TEveTrackPropagator *prop)
Set default propagator for tracks.
Definition: TEveTrack.cxx:627
void SelectByP(Double_t min_p, Double_t max_p)
Select visibility of tracks by momentum.
Definition: TEveTrack.cxx:1044
void SetRnrPoints(Bool_t r)
Set rendering of track as points for the list and the elements.
Definition: TEveTrack.cxx:772
Double_t fMaxPt
Definition: TEveTrack.h:156
void SanitizeMinMaxCuts()
Set Min/Max cuts so that they are within detected limits.
Definition: TEveTrack.cxx:727
Bool_t fRecurse
Definition: TEveTrack.h:150
Bool_t GetRnrLine() const
Definition: TEveTrack.h:199
virtual void SetLineStyle(Style_t s)
Set line style for the list and the elements.
Definition: TEveTrack.cxx:869
virtual void SetMarkerStyle(Style_t s)
Set marker style for the list and the elements.
Definition: TEveTrack.cxx:901
Double_t fMaxP
Definition: TEveTrack.h:159
Bool_t fRnrPoints
Definition: TEveTrack.h:153
TEveTrackPropagator * fPropagator
Definition: TEveTrack.h:148
TEveTrack * FindTrackByIndex(Int_t index)
Find track by index, select it and display it in the editor.
Definition: TEveTrack.cxx:1110
Bool_t GetRnrPoints() const
Definition: TEveTrack.h:203
virtual void SetLineColor(Color_t c)
Set the line color.
Definition: TEveTrack.h:183
TEveTrackList(const TEveTrackList &)
virtual void CopyVizParams(const TEveElement *el)
Copy visualization parameters from element el.
Definition: TEveTrack.cxx:1133
virtual void SetMarkerColor(Color_t c)
Set marker color for the list and the elements.
Definition: TEveTrack.cxx:933
void SelectByPt(Double_t min_pt, Double_t max_pt)
Select visibility of tracks by transverse momentum.
Definition: TEveTrack.cxx:999
virtual ~TEveTrackList()
Destructor.
Definition: TEveTrack.cxx:617
TEveTrack * FindTrackByLabel(Int_t label)
Find track by label, select it and display it in the editor.
Definition: TEveTrack.cxx:1087
void MakeTracks(Bool_t recurse=kTRUE)
Regenerate the visual representations of tracks.
Definition: TEveTrack.cxx:639
Double_t fLimP
Definition: TEveTrack.h:160
Double_t RoundMomentumLimit(Double_t x)
Round the momentum limit up to a nice value.
Definition: TEveTrack.cxx:714
virtual void SetMainColor(Color_t c)
Set main (line) color for the list and the elements.
Definition: TEveTrack.cxx:805
Double_t fLimPt
Definition: TEveTrack.h:157
Bool_t fRnrLine
Definition: TEveTrack.h:152
Double_t fMinP
Definition: TEveTrack.h:158
void SetRnrLine(Bool_t rnr)
Set rendering of track as line for the list and the elements.
Definition: TEveTrack.cxx:740
void FindMomentumLimits(TEveElement *el, Bool_t recurse=kTRUE)
Loop over track elements of argument el and find highest pT and p.
Definition: TEveTrack.cxx:696
virtual void SetMarkerSize(Size_t s)
Set marker size for the list and the elements.
Definition: TEveTrack.cxx:965
Holding structure for a number of track rendering parameters.
Bool_t GetFitCluster2Ds() const
virtual void GoToBounds(TEveVectorD &p)
Propagate particle to bounds.
void ResetTrack()
Reset cache holding particle trajectory.
Bool_t GetFitDaughters() const
Bool_t IntersectPlane(const TEveVectorD &p, const TEveVectorD &point, const TEveVectorD &normal, TEveVectorD &itsect)
Find intersection of currently propagated track with a plane.
Bool_t GetFitLineSegments() const
static Bool_t IsOutsideBounds(const TEveVectorD &point, Double_t maxRsqr, Double_t maxZ)
Double_t GetTrackLength(Int_t start_point=0, Int_t end_point=-1) const
Calculate track length from start_point to end_point.
Bool_t GetFitDecay() const
virtual Bool_t GoToVertex(TEveVectorD &v, TEveVectorD &p)
Propagate particle with momentum p to vertex v.
Int_t GetCurrentPoint() const
Get index of current point on track.
Double_t GetMaxZ() const
Bool_t GetFitReferences() const
void FillPointSet(TEvePointSet *ps) const
Reset ps and populate it with points in propagation cache.
Double_t GetMaxR() const
static TEveTrackPropagator fgDefault
virtual Bool_t GoToLineSegment(const TEveVectorD &s, const TEveVectorD &r, TEveVectorD &p)
Propagate particle with momentum p to line with start point s and vector r to the second point.
void InitTrack(const TEveVectorD &v, Int_t charge)
Initialize internal data-members for given particle parameters.
Visual representation of a track.
Definition: TEveTrack.h:33
vPathMark_t::iterator vPathMark_i
Definition: TEveTrack.h:43
virtual void WriteVizParams(std::ostream &out, const TString &var)
Write visualization parameters.
Definition: TEveTrack.cxx:491
virtual void SecSelected(TEveTrack *)
Emits "SecSelected(TEveTrack*)" signal.
Definition: TEveTrack.cxx:550
Int_t fCharge
Definition: TEveTrack.h:56
void SetPropagator(TEveTrackPropagator *prop)
Set track's render style.
Definition: TEveTrack.cxx:311
Int_t fPdg
Definition: TEveTrack.h:55
virtual ~TEveTrack()
Destructor.
Definition: TEveTrack.cxx:226
vPathMark_t & RefPathMarks()
Definition: TEveTrack.h:111
virtual void CopyVizParams(const TEveElement *el)
Copy visualization parameters from element el.
Definition: TEveTrack.cxx:477
virtual void MakeTrack(Bool_t recurse=kTRUE)
Calculate track representation based on track data and current settings of the propagator.
Definition: TEveTrack.cxx:340
TEveTrackPropagator * fPropagator
Last path-mark index tried in track-propagation.
Definition: TEveTrack.h:64
TEveVectorD fP
Definition: TEveTrack.h:51
virtual void ComputeBBox()
Compute the bounding box of the track.
Definition: TEveTrack.cxx:243
vPathMark_t::const_iterator vPathMark_ci
Definition: TEveTrack.h:44
void SetAttLineAttMarker(TEveTrackList *tl)
Set line and marker attributes from TEveTrackList.
Definition: TEveTrack.cxx:322
virtual TClass * ProjectedClass(const TEveProjection *p) const
Virtual from TEveProjectable, return TEveTrackProjected class.
Definition: TEveTrack.cxx:501
Int_t fIndex
Definition: TEveTrack.h:58
virtual void SetPathMarks(const TEveTrack &t)
Copy path-marks from t.
Definition: TEveTrack.cxx:301
void SortPathMarksByTime()
Sort registered pat-marks by time.
Definition: TEveTrack.cxx:518
void PrintPathMarks()
Print registered path-marks.
Definition: TEveTrack.cxx:526
virtual const TGPicture * GetListTreeIcon(Bool_t open=kFALSE)
Returns list-tree icon for TEveTrack.
Definition: TEveTrack.cxx:235
TEveTrack()
Default constructor.
Definition: TEveTrack.cxx:49
Double_t fDpDs
Definition: TEveTrack.h:54
Int_t fLabel
Definition: TEveTrack.h:57
Bool_t fLockPoints
Definition: TEveTrack.h:60
virtual void SetTrackParams(const TEveTrack &t)
Copy track parameters from t.
Definition: TEveTrack.cxx:284
Double_t fBeta
Definition: TEveTrack.h:53
TEveVectorD fV
Definition: TEveTrack.h:50
virtual void SetStdTitle()
Set standard track title based on most data-member values.
Definition: TEveTrack.cxx:268
vPathMark_t fPathMarks
Definition: TEveTrack.h:61
Int_t fLastPMIdx
Definition: TEveTrack.h:62
TEveVectorD fPEnd
Definition: TEveTrack.h:52
TT Perp2() const
Definition: TEveVector.h:100
TT Perp() const
Definition: TEveVector.h:101
TT Mag() const
Definition: TEveVector.h:98
virtual void * GetUserData() const =0
TGListTreeItem * GetSelected() const
Definition: TGListTree.h:397
void SetSelected(TGListTreeItem *item)
Definition: TGListTree.h:368
void HighlightItem(TGListTreeItem *item)
Highlight item.
virtual const char * GetName() const
Returns name of object.
Definition: TObject.cxx:357
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:866
Description of the static properties of a particle.
Definition: TParticlePDG.h:19
Int_t PdgCode() const
Definition: TParticlePDG.h:66
Double_t Charge() const
Definition: TParticlePDG.h:68
Description of the dynamic properties of a particle.
Definition: TParticle.h:26
TParticlePDG * GetPDG(Int_t mode=0) const
Returns a pointer to the TParticlePDG object using the pdgcode.
Definition: TParticle.cxx:273
virtual const char * GetName() const
Return particle name.
Definition: TParticle.cxx:257
virtual Int_t Size() const
Definition: TPolyMarker3D.h:73
virtual const char * GetName() const
Returns name of object.
Definition: TPolyMarker3D.h:57
virtual void SetName(const char *name)
Change (i.e.
Float_t * fP
Definition: TPolyMarker3D.h:36
void Emit(const char *signal, const T &arg)
Activate signal with single parameter.
Definition: TQObject.h:164
Basic string class.
Definition: TString.h:131
const char * Data() const
Definition: TString.h:364
Double_t x[n]
Definition: legend1.C:17
const Int_t n
Definition: legend1.C:16
static double P[]
TMath.
Definition: TMathBase.h:35
Int_t Nint(T x)
Round to nearest integer. Rounds half integers to the nearest even integer.
Definition: TMath.h:701
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:212
Double_t Floor(Double_t x)
Definition: TMath.h:691
Double_t Ceil(Double_t x)
Definition: TMath.h:683
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:723
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:180
Double_t Log10(Double_t x)
Definition: TMath.h:752
TCanvas * style()
Definition: style.C:1
auto * m
Definition: textangle.C:8
auto * a
Definition: textangle.C:12