CMS 3D CMS Logo

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