Logo ROOT  
Reference Guide
REveTrackPropagator.cxx
Go to the documentation of this file.
1// @(#)root/eve7:$Id$
2// Authors: Matevz Tadel & Alja Mrak-Tadel: 2006, 2007, 2018
3
4/*************************************************************************
5 * Copyright (C) 1995-2019, 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
13#include <ROOT/REveTrack.hxx>
14#include <ROOT/REveTrans.hxx>
15
16#include "TMath.h"
17
18using namespace ROOT::Experimental;
19namespace REX = ROOT::Experimental;
20
21/** \class REveMagField
22\ingroup REve
23Abstract base-class for interfacing to magnetic field needed by the
24REveTrackPropagator.
25See sub-classes for two simple implementations.
26
27NOTE: Magnetic field direction convention is inverted.
28*/
29
30
31/** \class REveMagFieldConst
32\ingroup REve
33Implements constant magnetic field, given by a vector fB.
34
35NOTE: Magnetic field direction convention is inverted.
36*/
37
38
39/** \class REveMagFieldDuo
40\ingroup REve
41Implements constant magnetic filed that switches on given axial radius fR2
42from vector fBIn to fBOut.
43
44NOTE: Magnetic field direction convention is inverted.
45*/
46
47namespace
48{
49 //const Double_t kBMin = 1e-6;
50 const Double_t kPtMinSqr = 1e-20;
51 const Double_t kAMin = 1e-10;
52 const Double_t kStepEps = 1e-3;
53}
54
55////////////////////////////////////////////////////////////////////////////////
56/// Default constructor.
57
58REveTrackPropagator::Helix_t::Helix_t() :
59 fCharge(0),
60 fMaxAng(45), fMaxStep(20.f), fDelta(0.1),
61 fPhi(0), fValid(kFALSE),
62 fLam(-1), fR(-1), fPhiStep(-1), fSin(-1), fCos(-1),
63 fRKStep(20.0),
64 fPtMag(-1), fPlMag(-1), fLStep(-1)
65{
66}
67
68////////////////////////////////////////////////////////////////////////////////
69/// Common update code for helix and RK propagation.
70
72{
73 fB = b;
74
75 // base vectors
76 fE1 = b;
77 fE1.Normalize();
78 fPlMag = p.Dot(fE1);
79 fPl = fE1*fPlMag;
80
81 fPt = p - fPl;
82 fPtMag = fPt.Mag();
83 fE2 = fPt;
84 fE2.Normalize();
85}
86
87////////////////////////////////////////////////////////////////////////////////
88/// Update helix parameters.
89
91 Bool_t full_update, Bool_t enforce_max_step)
92{
93 UpdateCommon(p, b);
94
95 // helix parameters
96 TMath::Cross(fE1.Arr(), fE2.Arr(), fE3.Arr());
97 if (fCharge < 0) fE3.NegateXYZ();
98
99 if (full_update)
100 {
101 using namespace TMath;
102 Double_t a = fgkB2C * b.Mag() * Abs(fCharge);
103 if (a > kAMin && fPtMag*fPtMag > kPtMinSqr)
104 {
105 fValid = kTRUE;
106
107 fR = Abs(fPtMag / a);
108 fLam = fPlMag / fPtMag;
109
110 // get phi step, compare fMaxAng with fDelta
111 fPhiStep = fMaxAng * DegToRad();
112 if (fR > fDelta)
113 {
114 Double_t ang = 2.0 * ACos(1.0f - fDelta/fR);
115 if (ang < fPhiStep)
116 fPhiStep = ang;
117 }
118
119 // check max step size
120 Double_t curr_step = fR * fPhiStep * Sqrt(1.0f + fLam*fLam);
121 if (curr_step > fMaxStep || enforce_max_step)
122 fPhiStep *= fMaxStep / curr_step;
123
124 fLStep = fR * fPhiStep * fLam;
125 fSin = Sin(fPhiStep);
126 fCos = Cos(fPhiStep);
127 }
128 else
129 {
130 fValid = kFALSE;
131 }
132 }
133}
134
135////////////////////////////////////////////////////////////////////////////////
136/// Update helix for stepper RungeKutta.
137
139{
140 UpdateCommon(p, b);
141
142 if (fCharge)
143 {
144 fValid = kTRUE;
145
146 // cached values for propagator
147 fB = b;
148 fPlMag = p.Dot(fB);
149 }
150 else
151 {
152 fValid = kFALSE;
153 }
154}
155
156////////////////////////////////////////////////////////////////////////////////
157/// Step helix for given momentum p from vertex v.
158
160 REveVector4D& vOut, REveVectorD& pOut)
161{
162 vOut = v;
163
164 if (fValid)
165 {
166 REveVectorD d = fE2*(fR*fSin) + fE3*(fR*(1-fCos)) + fE1*fLStep;
167 vOut += d;
168 vOut.fT += TMath::Abs(fLStep);
169
170 pOut = fPl + fE2*(fPtMag*fCos) + fE3*(fPtMag*fSin);
171
172 fPhi += fPhiStep;
173 }
174 else
175 {
176 // case: pT < kPtMinSqr or B < kBMin
177 // might happen if field directon changes pT ~ 0 or B becomes zero
178 vOut += p * (fMaxStep / p.Mag());
179 vOut.fT += fMaxStep;
180 pOut = p;
181 }
182}
183
184/** \class REveTrackPropagator
185\ingroup REve
186Holding structure for a number of track rendering parameters.
187Calculates path taking into account the parameters.
188
189NOTE: Magnetic field direction convention is inverted.
190
191This is decoupled from REveTrack/REveTrackList to allow sharing of the
192Propagator among several instances. Back references are kept so the tracks
193can be recreated when the parameters change.
194
195REveTrackList has Get/Set methods for RnrStlye.
196
197Enum EProjTrackBreaking_e and member fProjTrackBreaking specify whether 2D
198projected tracks get broken into several segments when the projected space
199consists of separate domains (like Rho-Z). The track-breaking is enabled by
200default.
201*/
202
204const Double_t REveTrackPropagator::fgkB2C = 0.299792458e-2;
206
209
210////////////////////////////////////////////////////////////////////////////////
211/// Default constructor.
212
213REveTrackPropagator::REveTrackPropagator(const std::string& n, const std::string& t,
214 REveMagField *field, Bool_t own_field) :
215 REveElement(n, t),
217
219 fMagFieldObj(field),
220 fOwnMagFiledObj(own_field),
221
222 fMaxR (350), fMaxZ (450),
223 fNMax (4096), fMaxOrbs (0.5),
224
231 fRnrFV (kFALSE),
232 fPMAtt(), fFVAtt(),
233
235
236 fV()
237{
241
245
249
250 if (!fMagFieldObj) {
253 }
254}
255
256////////////////////////////////////////////////////////////////////////////////
257/// Destructor.
258
260{
261 if (fOwnMagFiledObj)
262 {
263 delete fMagFieldObj;
264 }
265}
266
267////////////////////////////////////////////////////////////////////////////////
268/// Virtual from REveRefBackPtr - track reference count has reached zero.
269
271{
272 CheckReferenceCount("REveTrackPropagator::OnZeroRefCount ");
273}
274
275////////////////////////////////////////////////////////////////////////////////
276/// Check reference count - virtual from REveElement.
277/// Must also take into account references from REveRefBackPtr.
278
279void REveTrackPropagator::CheckReferenceCount(const std::string& from)
280{
281 if (fRefCount <= 0)
282 {
284 }
285}
286
287////////////////////////////////////////////////////////////////////////////////
288/// Element-change notification.
289/// Stamp all tracks as requiring display-list regeneration.
290
292{
293 for (auto &i: fBackRefs) {
294 auto track = dynamic_cast<REveTrack *>(i.first);
295 if (track) track->StampObjProps();
296 }
297}
298
299////////////////////////////////////////////////////////////////////////////////
300/// Initialize internal data-members for given particle parameters.
301
303{
304 fV = v;
305 fPoints.push_back(fV);
306
307 // init helix
308 fH.fPhi = 0;
309 fH.fCharge = charge;
310}
311
312////////////////////////////////////////////////////////////////////////////////
313/// REveVectorF wrapper.
314
316{
317 REveVectorD vd(v);
318 InitTrack(vd, charge);
319}
320
321////////////////////////////////////////////////////////////////////////////////
322/// Reset cache holding particle trajectory.
323
325{
326 fLastPoints.clear();
327 fPoints.swap(fLastPoints);
328
329 // reset helix
330 fH.fPhi = 0;
331}
332
333////////////////////////////////////////////////////////////////////////////////
334/// Get index of current point on track.
335
337{
338 return fPoints.size() - 1;
339}
340
341////////////////////////////////////////////////////////////////////////////////
342/// Calculate track length from start_point to end_point.
343/// If end_point is less than 0, distance to the end is returned.
344
346{
347 if (end_point < 0) end_point = fPoints.size() - 1;
348
349 Double_t sum = 0;
350 for (Int_t i = start_point; i < end_point; ++i)
351 {
352 sum += (fPoints[i+1] - fPoints[i]).Mag();
353 }
354 return sum;
355}
356
357////////////////////////////////////////////////////////////////////////////////
358/// Propagate particle with momentum p to vertex v.
359
361{
362 Update(fV, p, kTRUE);
363
364 if ((v-fV).Mag() < kStepEps)
365 {
366 fPoints.push_back(v);
367 return kTRUE;
368 }
369
370 return fH.fValid ? LoopToVertex(v, p) : LineToVertex(v);
371}
372
373////////////////////////////////////////////////////////////////////////////////
374/// Propagate particle with momentum p to line with start point s and vector r
375/// to the second point.
376
378{
379 Update(fV, p, kTRUE);
380
381 if (!fH.fValid)
382 {
386 return kTRUE;
387 }
388 else
389 {
390 return LoopToLineSegment(s, r, p);
391 }
392}
393
394////////////////////////////////////////////////////////////////////////////////
395/// REveVectorF wrapper.
396
398{
399 REveVectorD vd(v), pd(p);
400 Bool_t result = GoToVertex(vd, pd);
401 v = vd; p = pd;
402 return result;
403}
404
405////////////////////////////////////////////////////////////////////////////////
406/// REveVectorF wrapper.
407
409{
410 REveVectorD sd(s), rd(r), pd(p);
411 Bool_t result = GoToLineSegment(sd, rd, pd);
412 p = pd;
413 return result;
414}
415
416////////////////////////////////////////////////////////////////////////////////
417/// Propagate particle to bounds.
418/// Return TRUE if hit bounds.
419
421{
422 Update(fV, p, kTRUE);
423
425}
426
427////////////////////////////////////////////////////////////////////////////////
428/// REveVectorF wrapper.
429
431{
432 REveVectorD pd(p);
433 GoToBounds(pd);
434 p = pd;
435}
436
437////////////////////////////////////////////////////////////////////////////////
438/// Update helix / B-field projection state.
439
441 Bool_t full_update, Bool_t enforce_max_step)
442{
443 if (fStepper == kHelix)
444 {
445 fH.UpdateHelix(p, fMagFieldObj->GetFieldD(v), !fMagFieldObj->IsConst() || full_update, enforce_max_step);
446 }
447 else
448 {
450
451 if (full_update)
452 {
453 using namespace TMath;
454
456 if (a > kAMin)
457 {
458 fH.fR = p.Mag() / a;
459
460 // get phi step, compare fDelta with MaxAng
462 if (fH.fR > fH.fDelta )
463 {
464 Double_t ang = 2.0 * ACos(1.0f - fH.fDelta/fH.fR);
465 if (ang < fH.fPhiStep)
466 fH.fPhiStep = ang;
467 }
468
469 // check against maximum step-size
470 fH.fRKStep = fH.fR * fH.fPhiStep * Sqrt(1 + fH.fLam*fH.fLam);
471 if (fH.fRKStep > fH.fMaxStep || enforce_max_step)
472 {
475 }
476 }
477 else
478 {
480 }
481 }
482 }
483}
484
485////////////////////////////////////////////////////////////////////////////////
486/// Wrapper to step helix.
487
489{
490 if (fStepper == kHelix)
491 {
492 fH.Step(v, p, vOut, pOut);
493 }
494 else
495 {
496 Double_t vecRKIn[7];
497 vecRKIn[0] = v.fX;
498 vecRKIn[1] = v.fY;
499 vecRKIn[2] = v.fZ;
500 Double_t pm = p.Mag();
501 Double_t nm = 1.0 / pm;
502 vecRKIn[3] = p.fX*nm;
503 vecRKIn[4] = p.fY*nm;
504 vecRKIn[5] = p.fZ*nm;
505 vecRKIn[6] = p.Mag();
506
507 Double_t vecRKOut[7];
508 StepRungeKutta(fH.fRKStep, vecRKIn, vecRKOut);
509
510 vOut.fX = vecRKOut[0];
511 vOut.fY = vecRKOut[1];
512 vOut.fZ = vecRKOut[2];
513 vOut.fT = v.fT + fH.fRKStep;
514 pm = vecRKOut[6];
515 pOut.fX = vecRKOut[3]*pm;
516 pOut.fY = vecRKOut[4]*pm;
517 pOut.fZ = vecRKOut[5]*pm;
518 }
519}
520
521////////////////////////////////////////////////////////////////////////////////
522/// Propagate charged particle with momentum p to bounds.
523/// It is expected that Update() with full-update was called before.
524
526{
527 const Double_t maxRsq = fMaxR*fMaxR;
528
529 REveVector4D currV(fV);
530 REveVector4D forwV(fV);
531 REveVectorD forwP (p);
532
533 Int_t np = fPoints.size();
534 Double_t maxPhi = fMaxOrbs*TMath::TwoPi();
535
536 while (fH.fPhi < maxPhi && np<fNMax)
537 {
538 Step(currV, p, forwV, forwP);
539
540 // cross R
541 if (forwV.Perp2() > maxRsq)
542 {
543 Float_t t = (fMaxR - currV.R()) / (forwV.R() - currV.R());
544 if (t < 0 || t > 1)
545 {
546 Warning("HelixToBounds", "In MaxR crossing expected t>=0 && t<=1: t=%f, r1=%f, r2=%f, MaxR=%f.",
547 t, currV.R(), forwV.R(), fMaxR);
548 return;
549 }
550 REveVectorD d(forwV);
551 d -= currV;
552 d *= t;
553 d += currV;
554 fPoints.push_back(d);
555 return;
556 }
557
558 // cross Z
559 else if (TMath::Abs(forwV.fZ) > fMaxZ)
560 {
561 Double_t t = (fMaxZ - TMath::Abs(currV.fZ)) / TMath::Abs((forwV.fZ - currV.fZ));
562 if (t < 0 || t > 1)
563 {
564 Warning("HelixToBounds", "In MaxZ crossing expected t>=0 && t<=1: t=%f, z1=%f, z2=%f, MaxZ=%f.",
565 t, currV.fZ, forwV.fZ, fMaxZ);
566 return;
567 }
568 REveVectorD d(forwV -currV);
569 d *= t;
570 d += currV;
571 fPoints.push_back(d);
572 return;
573 }
574
575 currV = forwV;
576 p = forwP;
577 Update(currV, p);
578
579 fPoints.push_back(currV);
580 ++np;
581 }
582}
583
584////////////////////////////////////////////////////////////////////////////////
585/// Propagate charged particle with momentum p to vertex v.
586/// It is expected that Update() with full-update was called before.
587
589{
590 const Double_t maxRsq = fMaxR * fMaxR;
591
592 REveVector4D currV(fV);
593 REveVector4D forwV(fV);
594 REveVectorD forwP(p);
595
596 Int_t first_point = fPoints.size();
597 Int_t np = first_point;
598
599 Double_t prod0=0, prod1;
600
601 do
602 {
603 Step(currV, p, forwV, forwP);
604 Update(forwV, forwP);
605
606 if (PointOverVertex(v, forwV, &prod1))
607 {
608 break;
609 }
610
611 if (IsOutsideBounds(forwV, maxRsq, fMaxZ))
612 {
613 fV = currV;
614 return kFALSE;
615 }
616
617 fPoints.push_back(forwV);
618 currV = forwV;
619 p = forwP;
620 prod0 = prod1;
621 ++np;
622 } while (np < fNMax);
623
624 // make the remaining fractional step
625 if (np > first_point)
626 {
627 if ((v - currV).Mag() > kStepEps)
628 {
629 Double_t step_frac = prod0 / (prod0 - prod1);
630 if (step_frac > 0)
631 {
632 // Step for fraction of previous step size.
633 // We pass 'enforce_max_step' flag to Update().
634 Float_t orig_max_step = fH.fMaxStep;
635 fH.fMaxStep = step_frac * (forwV - currV).Mag();
636 Update(currV, p, kTRUE, kTRUE);
637 Step(currV, p, forwV, forwP);
638 p = forwP;
639 currV = forwV;
640 fPoints.push_back(currV);
641 ++np;
642 fH.fMaxStep = orig_max_step;
643 }
644
645 // Distribute offset to desired crossing point over all segment.
646
647 REveVectorD off(v - currV);
648 off *= 1.0f / currV.fT;
649 DistributeOffset(off, first_point, np, p);
650 fV = v;
651 return kTRUE;
652 }
653 }
654
655 fPoints.push_back(v);
656 fV = v;
657 return kTRUE;
658}
659
660////////////////////////////////////////////////////////////////////////////////
661/// Propagate charged particle with momentum p to line segment with point s and
662/// vector r to the second point. It is expected that Update() with full-update
663/// was called before. Returns kFALSE if hits bounds.
664
666{
667 const Double_t maxRsq = fMaxR * fMaxR;
668 const Double_t rMagInv = 1./r.Mag();
669
670 REveVector4D currV(fV);
671 REveVector4D forwV(fV);
672 REveVectorD forwP(p);
673
674 Int_t first_point = fPoints.size();
675 Int_t np = first_point;
676
677 REveVectorD forwC;
678 REveVectorD currC;
679 do
680 {
681 Step(currV, p, forwV, forwP);
682 Update(forwV, forwP);
683
684 ClosestPointFromVertexToLineSegment(forwV, s, r, rMagInv, forwC);
685
686 // check forwV is over segment with orthogonal component of
687 // momentum to vector r
688 REveVectorD b = r; b.Normalize();
689 Double_t x = forwP.Dot(b);
690 REveVectorD pTPM = forwP - x*b;
691 if (pTPM.Dot(forwC - forwV) < 0)
692 {
693 break;
694 }
695
696 if (IsOutsideBounds(forwV, maxRsq, fMaxZ))
697 {
698 fV = currV;
699 return kFALSE;
700 }
701
702 fPoints.push_back(forwV);
703 currV = forwV;
704 p = forwP;
705 currC = forwC;
706 ++np;
707 } while (np < fNMax);
708
709 // Get closest point on segment relative to line with forw and currV points.
711 ClosestPointBetweenLines(s, r, currV, forwV - currV, v);
712
713 // make the remaining fractional step
714 if (np > first_point)
715 {
716 if ((v - currV).Mag() > kStepEps)
717 {
718 REveVector last_step = forwV - currV;
719 REveVector delta = v - currV;
720 Double_t step_frac = last_step.Dot(delta) / last_step.Mag2();
721 if (step_frac > 0)
722 {
723 // Step for fraction of previous step size.
724 // We pass 'enforce_max_step' flag to Update().
725 Float_t orig_max_step = fH.fMaxStep;
726 fH.fMaxStep = step_frac * (forwV - currV).Mag();
727 Update(currV, p, kTRUE, kTRUE);
728 Step(currV, p, forwV, forwP);
729 p = forwP;
730 currV = forwV;
731 fPoints.push_back(currV);
732 ++np;
733 fH.fMaxStep = orig_max_step;
734 }
735
736 // Distribute offset to desired crossing point over all segment.
737
738 REveVectorD off(v - currV);
739 off *= 1.0f / currV.fT;
740 DistributeOffset(off, first_point, np, p);
741 fV = v;
742 return kTRUE;
743 }
744 }
745
746 fPoints.push_back(v);
747 fV = v;
748 return kTRUE;
749}
750
751////////////////////////////////////////////////////////////////////////////////
752/// Distribute offset between first and last point index and rotate
753/// momentum.
754
756{
757 // Calculate the required momentum rotation.
758 // lpd - last-points-delta
759 REveVectorD lpd0(fPoints[np-1]);
760 lpd0 -= fPoints[np-2];
761 lpd0.Normalize();
762
763 for (Int_t i = first_point; i < np; ++i)
764 {
765 fPoints[i] += off * fPoints[i].fT;
766 }
767
768 REveVectorD lpd1(fPoints[np-1]);
769 lpd1 -= fPoints[np-2];
770 lpd1.Normalize();
771
773 tt.SetupFromToVec(lpd0, lpd1);
774
775 // REveVectorD pb4(p);
776 // printf("Rotating momentum: p0 = "); p.Dump();
777 tt.RotateIP(p);
778 // printf(" p1 = "); p.Dump();
779 // printf(" n1=%f, n2=%f, dp = %f deg\n", pb4.Mag(), p.Mag(),
780 // TMath::RadToDeg()*TMath::ACos(p.Dot(pb4)/(pb4.Mag()*p.Mag())));
781}
782
783////////////////////////////////////////////////////////////////////////////////
784/// Propagate neutral particle to vertex v.
785
787{
788 REveVector4D currV = v;
789
790 currV.fX = v.fX;
791 currV.fY = v.fY;
792 currV.fZ = v.fZ;
793 fPoints.push_back(currV);
794
795 fV = v;
796 return kTRUE;
797}
798
799////////////////////////////////////////////////////////////////////////////////
800/// Propagate neutral particle with momentum p to bounds.
801
803{
804 Double_t tZ = 0, tR = 0, tB = 0;
805
806 // time where particle intersect +/- fMaxZ
807 if (p.fZ > 0)
808 tZ = (fMaxZ - fV.fZ) / p.fZ;
809 else if (p.fZ < 0)
810 tZ = - (fMaxZ + fV.fZ) / p.fZ;
811 else
812 tZ = 1e99;
813
814 // time where particle intersects cylinder
815 Double_t a = p.fX*p.fX + p.fY*p.fY;
816 Double_t b = 2.0 * (fV.fX*p.fX + fV.fY*p.fY);
818 Double_t d = b*b - 4.0*a*c;
819 if (d >= 0) {
820 Double_t sqrtD = TMath::Sqrt(d);
821 tR = (-b - sqrtD) / (2.0 * a);
822 if (tR < 0) {
823 tR = (-b + sqrtD) / (2.0 * a);
824 }
825 tB = tR < tZ ? tR : tZ; // compare the two times
826 } else {
827 tB = tZ;
828 }
829 REveVectorD nv(fV.fX + p.fX*tB, fV.fY + p.fY*tB, fV.fZ + p.fZ*tB);
830 LineToVertex(nv);
831}
832
833////////////////////////////////////////////////////////////////////////////////
834/// Intersect helix with a plane. Current position and argument p define
835/// the helix.
836
838 const REveVectorD& point,
839 const REveVectorD& normal,
840 REveVectorD& itsect)
841{
842 REveVectorD pos(fV);
843 REveVectorD mom(p);
844 if (fMagFieldObj->IsConst())
846
847 REveVectorD n(normal);
848 REveVectorD delta = pos - point;
849 Double_t d = delta.Dot(n);
850 if (d > 0) {
851 n.NegateXYZ(); // Turn normal around so that we approach from negative side of the plane
852 d = -d;
853 }
854
855 REveVector4D forwV;
856 REveVectorD forwP;
857 REveVector4D pos4(pos);
858 while (kTRUE)
859 {
860 Update(pos4, mom);
861 Step(pos4, mom, forwV , forwP);
862 Double_t new_d = (forwV - point).Dot(n);
863 if (new_d < d)
864 {
865 // We are going further away ... fail intersect.
866 Warning("HelixIntersectPlane", "going away from the plane.");
867 return kFALSE;
868 }
869 if (new_d > 0)
870 {
871 delta = forwV - pos;
872 itsect = pos + delta * (d / (d - new_d));
873 return kTRUE;
874 }
875 pos4 = forwV;
876 mom = forwP;
877 }
878}
879
880////////////////////////////////////////////////////////////////////////////////
881/// Intersect line with a plane. Current position and argument p define
882/// the line.
883
885 const REveVectorD& point,
886 const REveVectorD& normal,
887 REveVectorD& itsect)
888{
889 REveVectorD pos(fV.fX, fV.fY, fV.fZ);
890 REveVectorD delta = point - pos;
891
892 Double_t pn = p.Dot(normal);
893 if (pn == 0)
894 {
895 return kFALSE;
896 }
897 Double_t t = delta.Dot(normal) / pn;
898 if (t < 0) {
899 return kFALSE;
900 } else {
901 itsect = pos + p*t;
902 return kTRUE;
903 }
904}
905
906////////////////////////////////////////////////////////////////////////////////
907/// Find intersection of currently propagated track with a plane.
908/// Current track position is used as starting point.
909///
910/// Args:
911/// - p - track momentum to use for extrapolation
912/// - point - a point on a plane
913/// - normal - normal of the plane
914/// - itsect - output, point of intersection
915/// Returns:
916/// - kFALSE if intersection can not be found, kTRUE otherwise.
917
919 const REveVectorD& point,
920 const REveVectorD& normal,
921 REveVectorD& itsect)
922{
923 if (fH.fCharge && fMagFieldObj && p.Perp2() > kPtMinSqr)
924 return HelixIntersectPlane(p, point, normal, itsect);
925 else
926 return LineIntersectPlane(p, point, normal, itsect);
927}
928
929////////////////////////////////////////////////////////////////////////////////
930/// Get closest point from given vertex v to line segment defined with s and r.
931/// Argument rMagInv is cached. rMagInv= 1./rMag()
932
934 const REveVectorD& s,
935 const REveVectorD& r,
936 Double_t rMagInv,
937 REveVectorD& c)
938{
939 REveVectorD dir = v - s;
940 REveVectorD b1 = r * rMagInv;
941
942 // parallel distance
943 Double_t dot = dir.Dot(b1);
944 REveVectorD dirI = dot * b1;
945
946 Double_t facX = dot * rMagInv;
947
948 if (facX <= 0)
949 c = s;
950 else if (facX >= 1)
951 c = s + r;
952 else
953 c = s + dirI;
954}
955
956////////////////////////////////////////////////////////////////////////////////
957/// Get closest point on line defined with vector p0 and u.
958/// Return false if the point is forced on the line segment.
959
961 const REveVectorD& u,
962 const REveVectorD& q0,
963 const REveVectorD& v,
964 REveVectorD& out)
965{
966 REveVectorD w0 = p0 -q0;
967 Double_t a = u.Mag2();
968 Double_t b = u.Dot(v);
969 Double_t c = v.Mag2();
970 Double_t d = u.Dot(w0);
971 Double_t e = v.Dot(w0);
972
973 Double_t x = (b*e - c*d)/(a*c -b*b);
974 Bool_t force = (x < 0 || x > 1);
975 out = p0 + TMath::Range(0., 1., x) * u;
976 return force;
977}
978
979////////////////////////////////////////////////////////////////////////////////
980/// Reset ps and populate it with points in propagation cache.
981
983{
984 Int_t size = TMath::Min(fNMax, (Int_t)fPoints.size());
985 ps->Reset(size);
986 for (Int_t i = 0; i < size; ++i)
987 {
988 const REveVector4D& v = fPoints[i];
989 ps->SetNextPoint(v.fX, v.fY, v.fZ);
990 }
991}
992
993////////////////////////////////////////////////////////////////////////////////
994/// Rebuild all tracks using this render-style.
995
997{
998 for (auto &i: fBackRefs) {
999 auto track = dynamic_cast<REveTrack *>(i.first);
1000 if (track) {
1001 track->MakeTrack();
1002 track->StampObjProps();
1003 }
1004 }
1005}
1006
1007////////////////////////////////////////////////////////////////////////////////
1008/// Set constant magnetic field and rebuild tracks.
1009
1011{
1012 SetMagFieldObj(new REveMagFieldConst(bX, bY, bZ));
1013}
1014
1015////////////////////////////////////////////////////////////////////////////////
1016/// Set constant magnetic field and rebuild tracks.
1017
1019{
1021
1022 fMagFieldObj = field;
1023 fOwnMagFiledObj = own_field;
1024
1025 RebuildTracks();
1026}
1027
1028////////////////////////////////////////////////////////////////////////////////
1029
1031{
1033}
1034
1035////////////////////////////////////////////////////////////////////////////////
1036/// Set maximum radius and rebuild tracks.
1037
1039{
1040 fMaxR = x;
1041 RebuildTracks();
1042}
1043
1044////////////////////////////////////////////////////////////////////////////////
1045/// Set maximum z and rebuild tracks.
1046
1048{
1049 fMaxZ = x;
1050 RebuildTracks();
1051}
1052
1053////////////////////////////////////////////////////////////////////////////////
1054/// Set maximum number of orbits and rebuild tracks.
1055
1057{
1058 fMaxOrbs = x;
1059 RebuildTracks();
1060}
1061
1062////////////////////////////////////////////////////////////////////////////////
1063/// Set maximum step angle and rebuild tracks.
1064/// WARNING -- this method / variable was mis-named.
1065
1067{
1068 Warning("SetMinAng", "This method was mis-named, use SetMaxAng() instead!");
1069 SetMaxAng(x);
1070}
1071////////////////////////////////////////////////////////////////////////////////
1072/// Get maximum step angle.
1073/// WARNING -- this method / variable was mis-named.
1074
1076{
1077 Warning("GetMinAng", "This method was mis-named, use GetMaxAng() instead!");
1078 return GetMaxAng();
1079}
1080
1081////////////////////////////////////////////////////////////////////////////////
1082/// Set maximum step angle and rebuild tracks.
1083
1085{
1086 fH.fMaxAng = x;
1087 RebuildTracks();
1088}
1089
1090////////////////////////////////////////////////////////////////////////////////
1091/// Set maximum step-size and rebuild tracks.
1092
1094{
1095 fH.fMaxStep = x;
1096 RebuildTracks();
1097}
1098
1099////////////////////////////////////////////////////////////////////////////////
1100/// Set maximum error and rebuild tracks.
1101
1103{
1104 fH.fDelta = x;
1105 RebuildTracks();
1106}
1107
1108////////////////////////////////////////////////////////////////////////////////
1109/// Set daughter creation point fitting and rebuild tracks.
1110
1112{
1113 fFitDaughters = x;
1114 RebuildTracks();
1115}
1116
1117////////////////////////////////////////////////////////////////////////////////
1118/// Set track-reference fitting and rebuild tracks.
1119
1121{
1122 fFitReferences = x;
1123 RebuildTracks();
1124}
1125
1126////////////////////////////////////////////////////////////////////////////////
1127/// Set decay fitting and rebuild tracks.
1128
1130{
1131 fFitDecay = x;
1132 RebuildTracks();
1133}
1134
1135////////////////////////////////////////////////////////////////////////////////
1136/// Set line segment fitting and rebuild tracks.
1137
1139{
1141 RebuildTracks();
1142}
1143
1144////////////////////////////////////////////////////////////////////////////////
1145/// Set 2D-cluster fitting and rebuild tracks.
1146
1148{
1149 fFitCluster2Ds = x;
1150 RebuildTracks();
1151}
1152
1153////////////////////////////////////////////////////////////////////////////////
1154/// Set decay rendering and rebuild tracks.
1155
1157{
1158 fRnrDecay = rnr;
1159 RebuildTracks();
1160}
1161
1162////////////////////////////////////////////////////////////////////////////////
1163/// Set rendering of 2D-clusters and rebuild tracks.
1164
1166{
1167 fRnrCluster2Ds = rnr;
1168 RebuildTracks();
1169}
1170
1171////////////////////////////////////////////////////////////////////////////////
1172/// Set daughter rendering and rebuild tracks.
1173
1175{
1176 fRnrDaughters = rnr;
1177 RebuildTracks();
1178}
1179
1180////////////////////////////////////////////////////////////////////////////////
1181/// Set track-reference rendering and rebuild tracks.
1182
1184{
1185 fRnrReferences = rnr;
1186 RebuildTracks();
1187}
1188
1189////////////////////////////////////////////////////////////////////////////////
1190/// Set first-vertex rendering and rebuild tracks.
1191
1193{
1194 fRnrFV = x;
1195 RebuildTracks();
1196}
1197
1198////////////////////////////////////////////////////////////////////////////////
1199/// Set projection break-point mode and rebuild tracks.
1200
1202{
1204 RebuildTracks();
1205}
1206
1207////////////////////////////////////////////////////////////////////////////////
1208/// Set projection break-point rendering and rebuild tracks.
1209
1211{
1212 fRnrPTBMarkers = x;
1213 RebuildTracks();
1214}
1215
1216////////////////////////////////////////////////////////////////////////////////
1217/// Wrapper to step with method RungeKutta.
1218
1220 Double_t* vect, Double_t* vout)
1221{
1222 /// ******************************************************************
1223 /// * *
1224 /// * Runge-Kutta method for tracking a particle through a magnetic *
1225 /// * field. Uses Nystroem algorithm (See Handbook Nat. Bur. of *
1226 /// * Standards, procedure 25.5.20) *
1227 /// * *
1228 /// * Input parameters *
1229 /// * CHARGE Particle charge *
1230 /// * STEP Step size *
1231 /// * VECT Initial co-ords,direction cosines,momentum *
1232 /// * Output parameters *
1233 /// * VOUT Output co-ords,direction cosines,momentum *
1234 /// * User routine called *
1235 /// * CALL GUFLD(X,F) *
1236 /// * *
1237 /// * ==>Called by : <USER>, GUSWIM *
1238 /// * Authors R.Brun, M.Hansroul ********* *
1239 /// * V.Perevoztchikov (CUT STEP implementation) *
1240 /// * *
1241 /// * *
1242 /// ******************************************************************
1243
1244 Double_t h2, h4, f[4];
1245 Double_t /* xyzt[3], */ a, b, c, ph,ph2;
1246 Double_t secxs[4],secys[4],seczs[4],hxp[3];
1247 Double_t g1, g2, g3, g4, g5, g6, ang2, dxt, dyt, dzt;
1248 Double_t est, at, bt, ct, cba;
1249 Double_t f1, f2, f3, f4, rho, tet, hnorm, hp, rho1, sint, cost;
1250
1251 Double_t x;
1252 Double_t y;
1253 Double_t z;
1254
1255 Double_t xt;
1256 Double_t yt;
1257 Double_t zt;
1258
1259 // const Int_t maxit = 1992;
1260 const Int_t maxit = 500;
1261 const Int_t maxcut = 11;
1262
1263 const Double_t hmin = 1e-4; // !!! MT ADD, should be member
1264 const Double_t kdlt = 1e-3; // !!! MT CHANGE from 1e-4, should be member
1265 const Double_t kdlt32 = kdlt/32.;
1266 const Double_t kthird = 1./3.;
1267 const Double_t khalf = 0.5;
1268 const Double_t kec = 2.9979251e-3;
1269
1270 const Double_t kpisqua = 9.86960440109;
1271 const Int_t kix = 0;
1272 const Int_t kiy = 1;
1273 const Int_t kiz = 2;
1274 const Int_t kipx = 3;
1275 const Int_t kipy = 4;
1276 const Int_t kipz = 5;
1277
1278 // *.
1279 // *. ------------------------------------------------------------------
1280 // *.
1281 // * this constant is for units cm,gev/c and kgauss
1282 // *
1283 Int_t iter = 0;
1284 Int_t ncut = 0;
1285 for(Int_t j = 0; j < 7; j++)
1286 vout[j] = vect[j];
1287
1288 Double_t pinv = kec * fH.fCharge / vect[6];
1289 Double_t tl = 0.;
1290 Double_t h = step;
1291 Double_t rest;
1292
1293 do {
1294 rest = step - tl;
1295 if (TMath::Abs(h) > TMath::Abs(rest))
1296 h = rest;
1297
1298 f[0] = -fH.fB.fX;
1299 f[1] = -fH.fB.fY;
1300 f[2] = -fH.fB.fZ;
1301
1302 // * start of integration
1303 x = vout[0];
1304 y = vout[1];
1305 z = vout[2];
1306 a = vout[3];
1307 b = vout[4];
1308 c = vout[5];
1309
1310 h2 = khalf * h;
1311 h4 = khalf * h2;
1312 ph = pinv * h;
1313 ph2 = khalf * ph;
1314 secxs[0] = (b * f[2] - c * f[1]) * ph2;
1315 secys[0] = (c * f[0] - a * f[2]) * ph2;
1316 seczs[0] = (a * f[1] - b * f[0]) * ph2;
1317 ang2 = (secxs[0]*secxs[0] + secys[0]*secys[0] + seczs[0]*seczs[0]);
1318 if (ang2 > kpisqua) break;
1319
1320 dxt = h2 * a + h4 * secxs[0];
1321 dyt = h2 * b + h4 * secys[0];
1322 dzt = h2 * c + h4 * seczs[0];
1323 xt = x + dxt;
1324 yt = y + dyt;
1325 zt = z + dzt;
1326
1327 // * second intermediate point
1328 est = TMath::Abs(dxt) + TMath::Abs(dyt) + TMath::Abs(dzt);
1329 if (est > h) {
1330 if (ncut++ > maxcut) break;
1331 h *= khalf;
1332 continue;
1333 }
1334
1335 // xyzt[0] = xt;
1336 // xyzt[1] = yt;
1337 // xyzt[2] = zt;
1338
1339 fH.fB = fMagFieldObj->GetFieldD(xt, yt, zt);
1340 f[0] = -fH.fB.fX;
1341 f[1] = -fH.fB.fY;
1342 f[2] = -fH.fB.fZ;
1343
1344 at = a + secxs[0];
1345 bt = b + secys[0];
1346 ct = c + seczs[0];
1347
1348 secxs[1] = (bt * f[2] - ct * f[1]) * ph2;
1349 secys[1] = (ct * f[0] - at * f[2]) * ph2;
1350 seczs[1] = (at * f[1] - bt * f[0]) * ph2;
1351 at = a + secxs[1];
1352 bt = b + secys[1];
1353 ct = c + seczs[1];
1354 secxs[2] = (bt * f[2] - ct * f[1]) * ph2;
1355 secys[2] = (ct * f[0] - at * f[2]) * ph2;
1356 seczs[2] = (at * f[1] - bt * f[0]) * ph2;
1357 dxt = h * (a + secxs[2]);
1358 dyt = h * (b + secys[2]);
1359 dzt = h * (c + seczs[2]);
1360 xt = x + dxt;
1361 yt = y + dyt;
1362 zt = z + dzt;
1363 at = a + 2.*secxs[2];
1364 bt = b + 2.*secys[2];
1365 ct = c + 2.*seczs[2];
1366
1367 est = TMath::Abs(dxt)+TMath::Abs(dyt)+TMath::Abs(dzt);
1368 if (est > 2.*TMath::Abs(h)) {
1369 if (ncut++ > maxcut) break;
1370 h *= khalf;
1371 continue;
1372 }
1373
1374 // xyzt[0] = xt;
1375 // xyzt[1] = yt;
1376 // xyzt[2] = zt;
1377
1378 fH.fB = fMagFieldObj->GetFieldD(xt, yt, zt);
1379 f[0] = -fH.fB.fX;
1380 f[1] = -fH.fB.fY;
1381 f[2] = -fH.fB.fZ;
1382
1383 z = z + (c + (seczs[0] + seczs[1] + seczs[2]) * kthird) * h;
1384 y = y + (b + (secys[0] + secys[1] + secys[2]) * kthird) * h;
1385 x = x + (a + (secxs[0] + secxs[1] + secxs[2]) * kthird) * h;
1386
1387 secxs[3] = (bt*f[2] - ct*f[1])* ph2;
1388 secys[3] = (ct*f[0] - at*f[2])* ph2;
1389 seczs[3] = (at*f[1] - bt*f[0])* ph2;
1390 a = a+(secxs[0]+secxs[3]+2. * (secxs[1]+secxs[2])) * kthird;
1391 b = b+(secys[0]+secys[3]+2. * (secys[1]+secys[2])) * kthird;
1392 c = c+(seczs[0]+seczs[3]+2. * (seczs[1]+seczs[2])) * kthird;
1393
1394 est = TMath::Abs(secxs[0]+secxs[3] - (secxs[1]+secxs[2]))
1395 + TMath::Abs(secys[0]+secys[3] - (secys[1]+secys[2]))
1396 + TMath::Abs(seczs[0]+seczs[3] - (seczs[1]+seczs[2]));
1397
1398 if (est > kdlt && TMath::Abs(h) > hmin) {
1399 if (ncut++ > maxcut) break;
1400 h *= khalf;
1401 continue;
1402 }
1403
1404 ncut = 0;
1405 // * if too many iterations, go to helix
1406 if (iter++ > maxit) break;
1407
1408 tl += h;
1409 if (est < kdlt32)
1410 h *= 2.;
1411 cba = 1./ TMath::Sqrt(a*a + b*b + c*c);
1412 vout[0] = x;
1413 vout[1] = y;
1414 vout[2] = z;
1415 vout[3] = cba*a;
1416 vout[4] = cba*b;
1417 vout[5] = cba*c;
1418 rest = step - tl;
1419 if (step < 0.) rest = -rest;
1420 if (rest < 1.e-5*TMath::Abs(step))
1421 {
1422 Float_t dot = (vout[3]*vect[3] + vout[4]*vect[4] + vout[5]*vect[5]);
1423 fH.fPhi += TMath::ACos(dot);
1424 return;
1425 }
1426
1427 } while(1);
1428
1429 // angle too big, use helix
1430
1431 f1 = f[0];
1432 f2 = f[1];
1433 f3 = f[2];
1434 f4 = TMath::Sqrt(f1*f1+f2*f2+f3*f3);
1435 rho = -f4*pinv;
1436 tet = rho * step;
1437
1438 hnorm = 1./f4;
1439 f1 = f1*hnorm;
1440 f2 = f2*hnorm;
1441 f3 = f3*hnorm;
1442
1443 hxp[0] = f2*vect[kipz] - f3*vect[kipy];
1444 hxp[1] = f3*vect[kipx] - f1*vect[kipz];
1445 hxp[2] = f1*vect[kipy] - f2*vect[kipx];
1446
1447 hp = f1*vect[kipx] + f2*vect[kipy] + f3*vect[kipz];
1448
1449 rho1 = 1./rho;
1450 sint = TMath::Sin(tet);
1451 cost = 2.*TMath::Sin(khalf*tet)*TMath::Sin(khalf*tet);
1452
1453 g1 = sint*rho1;
1454 g2 = cost*rho1;
1455 g3 = (tet-sint) * hp*rho1;
1456 g4 = -cost;
1457 g5 = sint;
1458 g6 = cost * hp;
1459
1460 vout[kix] = vect[kix] + g1*vect[kipx] + g2*hxp[0] + g3*f1;
1461 vout[kiy] = vect[kiy] + g1*vect[kipy] + g2*hxp[1] + g3*f2;
1462 vout[kiz] = vect[kiz] + g1*vect[kipz] + g2*hxp[2] + g3*f3;
1463
1464 vout[kipx] = vect[kipx] + g4*vect[kipx] + g5*hxp[0] + g6*f1;
1465 vout[kipy] = vect[kipy] + g4*vect[kipy] + g5*hxp[1] + g6*f2;
1466 vout[kipz] = vect[kipz] + g4*vect[kipz] + g5*hxp[2] + g6*f3;
1467
1468 fH.fPhi += tet;
1469}
ROOT::R::TRInterface & r
Definition: Object.C:4
#define d(i)
Definition: RSha256.hxx:102
#define b(i)
Definition: RSha256.hxx:100
#define f(i)
Definition: RSha256.hxx:104
#define c(i)
Definition: RSha256.hxx:101
#define h(i)
Definition: RSha256.hxx:106
#define e(i)
Definition: RSha256.hxx:103
int Int_t
Definition: RtypesCore.h:41
unsigned char UChar_t
Definition: RtypesCore.h:34
const Bool_t kFALSE
Definition: RtypesCore.h:88
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
float Float_t
Definition: RtypesCore.h:53
const Bool_t kTRUE
Definition: RtypesCore.h:87
@ kRed
Definition: Rtypes.h:64
@ kBlue
Definition: Rtypes.h:64
@ kYellow
Definition: Rtypes.h:64
void Warning(const char *location, const char *msgfmt,...)
virtual void CheckReferenceCount(const std::string &from="<unknown>")
Check external references to this and eventually auto-destruct the render-element.
REveMagFieldConst Interface to constant magnetic field.
REveMagField Abstract interface to magnetic field.
virtual Double_t GetMaxFieldMagD() const
REveVectorD GetFieldD(const REveVectorD &v) const
virtual void PrintField(Double_t x, Double_t y, Double_t z) const
REveRefBackPtr reference-count with back pointers.
Definition: REveUtil.hxx:131
REveTrackPropagator Calculates path of a particle taking into account special path-marks and imposed ...
virtual Bool_t GoToVertex(REveVectorD &v, REveVectorD &p)
Propagate particle with momentum p to vertex v.
REveTrackPropagator(const REveTrackPropagator &)
void Update(const REveVector4D &v, const REveVectorD &p, Bool_t full_update=kFALSE, Bool_t enforce_max_step=kFALSE)
Update helix / B-field projection state.
void SetRnrPTBMarkers(Bool_t x)
Set projection break-point rendering and rebuild tracks.
static Bool_t IsOutsideBounds(const REveVectorD &point, Double_t maxRsqr, Double_t maxZ)
void SetRnrDaughters(Bool_t x)
Set daughter rendering and rebuild tracks.
Bool_t HelixIntersectPlane(const REveVectorD &p, const REveVectorD &point, const REveVectorD &normal, REveVectorD &itsect)
Intersect helix with a plane.
void SetFitLineSegments(Bool_t x)
Set line segment fitting and rebuild tracks.
void SetRnrCluster2Ds(Bool_t x)
Set rendering of 2D-clusters and rebuild tracks.
Bool_t LineIntersectPlane(const REveVectorD &p, const REveVectorD &point, const REveVectorD &normal, REveVectorD &itsect)
Intersect line with a plane.
void CheckReferenceCount(const std::string &from="<unknown>") override
Check reference count - virtual from REveElement.
Double_t GetMinAng() const
Get maximum step angle.
void ResetTrack()
Reset cache holding particle trajectory.
void Step(const REveVector4D &v, const REveVectorD &p, REveVector4D &vOut, REveVectorD &pOut)
Wrapper to step helix.
Double_t GetTrackLength(Int_t start_point=0, Int_t end_point=-1) const
Calculate track length from start_point to end_point.
void InitTrack(const REveVectorD &v, Int_t charge)
Initialize internal data-members for given particle parameters.
void LineToBounds(REveVectorD &p)
Propagate neutral particle with momentum p to bounds.
void SetRnrReferences(Bool_t x)
Set track-reference rendering and rebuild tracks.
void ClosestPointFromVertexToLineSegment(const REveVectorD &v, const REveVectorD &s, const REveVectorD &r, Double_t rMagInv, REveVectorD &c)
Get closest point from given vertex v to line segment defined with s and r.
Bool_t IntersectPlane(const REveVectorD &p, const REveVectorD &point, const REveVectorD &normal, REveVectorD &itsect)
Find intersection of currently propagated track with a plane.
void RebuildTracks()
Rebuild all tracks using this render-style.
void SetMaxAng(Double_t x)
Set maximum step angle and rebuild tracks.
Int_t GetCurrentPoint() const
Get index of current point on track.
void SetMaxZ(Double_t x)
Set maximum z and rebuild tracks.
void SetMaxStep(Double_t x)
Set maximum step-size and rebuild tracks.
void OnZeroRefCount() override
Virtual from REveRefBackPtr - track reference count has reached zero.
void SetRnrDecay(Bool_t x)
Set decay rendering and rebuild tracks.
Bool_t PointOverVertex(const REveVector4D &v0, const REveVector4D &v, Double_t *p=0)
void FillPointSet(REvePointSet *ps) const
Reset ps and populate it with points in propagation cache.
Bool_t LoopToLineSegment(const REveVectorD &s, const REveVectorD &r, REveVectorD &p)
Propagate charged particle with momentum p to line segment with point s and vector r to the second po...
void DistributeOffset(const REveVectorD &off, Int_t first_point, Int_t np, REveVectorD &p)
Distribute offset between first and last point index and rotate momentum.
void SetFitDaughters(Bool_t x)
Set daughter creation point fitting and rebuild tracks.
void SetFitDecay(Bool_t x)
Set decay fitting and rebuild tracks.
void SetFitReferences(Bool_t x)
Set track-reference fitting and rebuild tracks.
void StampAllTracks()
Element-change notification.
void SetMaxOrbs(Double_t x)
Set maximum number of orbits and rebuild tracks.
void StepRungeKutta(Double_t step, Double_t *vect, Double_t *vout)
Wrapper to step with method RungeKutta.
void SetProjTrackBreaking(UChar_t x)
Set projection break-point mode and rebuild tracks.
Bool_t ClosestPointBetweenLines(const REveVectorD &, const REveVectorD &, const REveVectorD &, const REveVectorD &, REveVectorD &out)
Get closest point on line defined with vector p0 and u.
void SetMagFieldObj(REveMagField *field, Bool_t own_field=kTRUE)
Set constant magnetic field and rebuild tracks.
void SetRnrFV(Bool_t x)
Set first-vertex rendering and rebuild tracks.
void LoopToBounds(REveVectorD &p)
Propagate charged particle with momentum p to bounds.
void PrintMagField(Double_t x, Double_t y, Double_t z) const
void SetMaxR(Double_t x)
Set maximum radius and rebuild tracks.
virtual Bool_t GoToLineSegment(const REveVectorD &s, const REveVectorD &r, REveVectorD &p)
Propagate particle with momentum p to line with start point s and vector r to the second point.
void SetDelta(Double_t x)
Set maximum error and rebuild tracks.
void SetMinAng(Double_t x)
Set maximum step angle and rebuild tracks.
virtual void GoToBounds(REveVectorD &p)
Propagate particle to bounds.
Bool_t LoopToVertex(REveVectorD &v, REveVectorD &p)
Propagate charged particle with momentum p to vertex v.
void SetMagField(Double_t bX, Double_t bY, Double_t bZ)
Set constant magnetic field and rebuild tracks.
Bool_t LineToVertex(REveVectorD &v)
Propagate neutral particle to vertex v.
void SetFitCluster2Ds(Bool_t x)
Set 2D-cluster fitting and rebuild tracks.
REveTrack Track with given vertex, momentum and optional referece-points (path-marks) along its path.
Definition: REveTrack.hxx:40
virtual void MakeTrack(Bool_t recurse=kTRUE)
Calculate track representation based on track data and current settings of the propagator.
Definition: REveTrack.cxx:327
REveVector4T A four-vector template without TObject inheritance and virtual functions.
Definition: REveVector.hxx:239
TT Normalize(TT length=1)
Normalize the vector to length if current length is non-zero.
Definition: REveVector.cxx:56
TT Dot(const REveVectorT &a) const
Definition: REveVector.hxx:164
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition: TAttMarker.h:38
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition: TAttMarker.h:40
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition: TAttMarker.h:41
T Mag(const SVector< T, D > &rhs)
Vector magnitude (Euclidian norm) Compute : .
Definition: Functions.h:252
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
const Int_t n
Definition: legend1.C:16
TF1 * f1
Definition: legend1.C:11
Double_t Sqrt(Double_t x)
static constexpr double nm
static constexpr double s
static constexpr double ps
TMath.
Definition: TMathBase.h:35
Double_t ACos(Double_t)
Definition: TMath.h:658
Short_t Range(Short_t lb, Short_t ub, Short_t x)
Definition: TMathBase.h:244
constexpr Double_t DegToRad()
Conversion from degree to radian:
Definition: TMath.h:82
Double_t Sqrt(Double_t x)
Definition: TMath.h:681
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:180
Double_t Cos(Double_t)
Definition: TMath.h:631
Double_t Sin(Double_t)
Definition: TMath.h:627
T * Cross(const T v1[3], const T v2[3], T out[3])
Calculate the Cross Product of two vectors: out = [v1 x v2].
Definition: TMath.h:1165
Short_t Abs(Short_t d)
Definition: TMathBase.h:120
constexpr Double_t TwoPi()
Definition: TMath.h:45
#define Dot(u, v)
Definition: normal.c:49
void UpdateHelix(const REveVectorD &p, const REveVectorD &b, Bool_t full_update, Bool_t enforce_max_step)
Update helix parameters.
void Step(const REveVector4D &v, const REveVectorD &p, REveVector4D &vOut, REveVectorD &pOut)
Step helix for given momentum p from vertex v.
void UpdateRK(const REveVectorD &p, const REveVectorD &b)
Update helix for stepper RungeKutta.
void UpdateCommon(const REveVectorD &p, const REveVectorD &b)
Common update code for helix and RK propagation.
auto * tt
Definition: textangle.C:16
auto * a
Definition: textangle.C:12
static long int sum(long int i)
Definition: Factory.cxx:2276