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