CMS 3D CMS Logo

TwoTrackMinimumDistanceHelixHelix.cc
Go to the documentation of this file.
7 
8 using namespace std;
9 
11 theH(nullptr), theG(nullptr), pointsUpdated(false), themaxjump(20),thesingjacI(1./0.1), themaxiter(4)
12 { }
13 
15 
17  const GlobalPoint & gpH, const GlobalPoint & gpG ) {
18 
19  const double Bc2kH = theH->magneticField().inInverseGeV(gpH).z();
20  const double Bc2kG = theG->magneticField().inInverseGeV(gpG).z();
21  const double Ht= theH->momentum().perp();
22  const double Gt= theG->momentum().perp();
23  // thelambdaH=asin ( theH->momentum().z() / Hn );
24 
25  if ( Ht == 0. || Gt == 0. ) {
26  edm::LogWarning ("TwoTrackMinimumDistanceHelixHelix")
27  << "transverse momentum of input trajectory is zero.";
28  return true;
29  };
30 
31  if ( theH->charge() == 0. || theG->charge() == 0. ) {
32  edm::LogWarning ("TwoTrackMinimumDistanceHelixHelix")
33  << "charge of input track is zero.";
34  return true;
35  };
36 
37  if ( Bc2kG == 0.) {
38  edm::LogWarning ("TwoTrackMinimumDistanceHelixHelix")
39  << "magnetic field at point " << gpG << " is zero.";
40  return true;
41  };
42 
43  if ( Bc2kH == 0. ) {
44  edm::LogWarning ("TwoTrackMinimumDistanceHelixHelix")
45  << "magnetic field at point " << gpH << " is zero.";
46  return true;
47  };
48 
49  theh= Ht / (theH->charge() * Bc2kH );
50  thesinpH0= - theH->momentum().y() / ( Ht );
51  thecospH0= - theH->momentum().x() / ( Ht );
52  thetanlambdaH = - theH->momentum().z() / ( Ht);
53  thepH0 = atan2 ( thesinpH0 , thecospH0 );
54 
55  // thelambdaG=asin ( theG->momentum().z()/( Gn ) );
56 
57  theg= Gt / (theG->charge() * Bc2kG );
58  thesinpG0= - theG->momentum().y() / ( Gt);
59  thecospG0= - theG->momentum().x() / (Gt);
60  thetanlambdaG = - theG->momentum().z() / ( Gt);
61  thepG0 = atan2 ( thesinpG0 , thecospG0 );
62 
63  thea=theH->position().x() - theG->position().x() + theg * thesinpG0 -
64  theh * thesinpH0;
65  theb=theH->position().y() - theG->position().y() - theg * thecospG0 +
66  theh * thecospH0;
69  thed1= -theg * thetanlambdaG * thetanlambdaH;
70  thed2= theh * thetanlambdaG * thetanlambdaH;
71  thee1= thetanlambdaH * ( theH->position().z() - theG->position().z() -
72  theh * thepH0 * thetanlambdaH + theg * thepG0 * thetanlambdaG );
73  thee2= thetanlambdaG * ( theH->position().z() - theG->position().z() -
74  theh * thepH0 * thetanlambdaH + theg * thepG0 * thetanlambdaG );
75  return false;
76 }
77 
78 bool TwoTrackMinimumDistanceHelixHelix::oneIteration(double & dH, double & dG ) {
83 
84  const double A11= theh * ( - thesinpH * ( thea - theg * thesinpG ) +
85  thecospH * ( theb + theg * thecospG ) + thec1);
86  if (A11 < 0) { return true; };
87  const double A22= -theg * (- thesinpG * ( thea + theh * thesinpH ) +
88  thecospG*( theb - theh * thecospH ) + thec2);
89  if (A22 < 0) { return true; };
90  const double A12= theh * (-theg * thecospG * thecospH -
91  theg * thesinpH * thesinpG + thed1);
92  const double A21= -theg * (theh * thecospG * thecospH +
93  theh * thesinpH * thesinpG + thed2);
94  const double detaI = 1./(A11 * A22 - A12 * A21);
95  const double z1=theh * ( thecospH * ( thea - theg * thesinpG ) + thesinpH *
96  ( theb + theg*thecospG ) + thec1 * thepH + thed1 * thepG + thee1);
97  const double z2=-theg * (thecospG * ( thea + theh * thesinpH ) + thesinpG *
98  ( theb - theh*thecospH ) + thec2 * thepG + thed2 * thepH + thee2);
99 
100  dH=( z1 * A22 - z2 * A12 ) * detaI;
101  dG=( z2 * A11 - z1 * A21 ) * detaI;
102  if ( fabs(detaI) > thesingjacI ) { return true; };
103 
104  thepH -= dH;
105  thepG -= dG;
106  return false;
107 }
108 
109 /*
110 bool TwoTrackMinimumDistanceHelixHelix::parallelTracks() const {
111  return (fabs(theH->momentum().x() - theG->momentum().x()) < .00000001 )
112  && (fabs(theH->momentum().y() - theG->momentum().y()) < .00000001 )
113  && (fabs(theH->momentum().z() - theG->momentum().z()) < .00000001 )
114  && (theH->charge()==theG->charge())
115  ;
116 }
117 */
118 
121  const GlobalTrajectoryParameters & H, const float qual ) {
122  pointsUpdated = false;
123  theG= &G;
124  theH= &H;
125  bool retval=false;
126 
127  if ( updateCoeffs ( theG->position(), theH->position() ) ){
128  finalPoints();
129  return true;
130  }
131 
132  thepG = thepG0;
133  thepH = thepH0;
134 
135  int counter=0;
136  double pH=0; double pG=0;
137  do {
138  retval=oneIteration ( pG, pH );
139  if ( edm::isNotFinite(pG) || edm::isNotFinite(pH) ) retval=true;
140  if ( counter++>themaxiter ) retval=true;
141  } while ( (!retval) && ( fabs(pG) > qual || fabs(pH) > qual ));
142  if ( fabs ( theg * ( thepG - thepG0 ) ) > themaxjump ) retval=true;
143  if ( fabs ( theh * ( thepH - thepH0 ) ) > themaxjump ) retval=true;
144 
145  finalPoints();
146 
147  return retval;
148 }
149 
150 
152  if (pointsUpdated) return;
153  GlobalVector tmpG( sin(thepG) - thesinpG0,
154  - cos(thepG) + thecospG0,
155  thetanlambdaG * ( thepG- thepG0 )
156  );
157  pointG = theG->position() + theg * tmpG;
159 
160  GlobalVector tmpH( sin(thepH) - thesinpH0,
161  - cos(thepH) + thecospH0,
162  thetanlambdaH * ( thepH- thepH0 )
163  );
164  pointH = theH->position() + theh * tmpH;
166 
167  pointsUpdated = true;
168 }
169 
T perp() const
Definition: PV3DBase.h:72
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
#define nullptr
T y() const
Definition: PV3DBase.h:63
GlobalVector inInverseGeV(const GlobalPoint &gp) const
Field value ad specified global point, in 1/Gev.
Definition: MagneticField.h:41
bool updateCoeffs(const GlobalPoint &, const GlobalPoint &)
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
bool calculate(const GlobalTrajectoryParameters &, const GlobalTrajectoryParameters &, const float qual=.001)
const MagneticField & magneticField() const
T x() const
Definition: PV3DBase.h:62