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