CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ThirdHitPredictionFromCircle.cc
Go to the documentation of this file.
1 #include <cmath>
4 
7 
9 
10 // there are tons of safety checks.
11 // Try to move all of the out the regular control flow using gcc magic
12 
13 #if defined(__GNUC__) && (__GNUC__ > 2) && defined(__OPTIMIZE__)
14 # define likely(expr) (__builtin_expect(!!(expr), 1))
15 # define unlikely(expr) (__builtin_expect(!!(expr), 0))
16 #else
17 # define likely(expr) (expr)
18 # define unlikely(expr) (expr)
19 #endif
20 
23 
24 namespace {
25  template<class T> static inline T sqr(T t) { return t * t; }
26  template<class T> static inline T sgn(T t) { return std::signbit(t) ? -1. : 1.; }
27  template<class T> static inline T clamped_acos(T t)
28  { return unlikely(t <= -1) ? M_PI : unlikely(t >= 1) ? T(0) : std::acos(t); }
29  template<class T> static inline T clamped_sqrt(T t)
30  { return likely(t > 0) ? std::sqrt(t) : T(0); }
31 }
32 
34  const GlobalPoint& P1, const GlobalPoint& P2, float tolerance)
35  : p1(P1.x(), P1.y()), theTolerance(tolerance)
36 {
37  Point2D p2(P2.x(), P2.y());
38  Point2D diff = 0.5 * (p2 - p1);
39  delta2 = diff.mag2();
41  axis = Point2D(-diff.y(), diff.x()) / delta;
42  center = p1 + diff;
43 }
44 
46 {
47  double phi;
48  if (unlikely(std::abs(curvature) < 1.0e-5)) {
49  double cos = (center * axis) / radius;
50  phi = axis.phi() - clamped_acos(cos);
51  } else {
52  double sign = sgn(curvature);
53  double radius2 = sqr(1.0 / curvature);
54  double orthog = clamped_sqrt(radius2 - delta2);
55  Basic2DVector<double> lcenter = center - sign * orthog * axis;
56  double rc2 = lcenter.mag2();
57  double cos = (rc2 + sqr(radius) - radius2) /
58  (2. *std:: sqrt(rc2) * radius);
59  phi = lcenter.phi() + sign * clamped_acos(cos);
60  }
61 
62  while(unlikely(phi >= M_PI)) phi -= 2. * M_PI;
63  while(unlikely(phi < -M_PI)) phi += 2. * M_PI;
64 
65  return phi;
66 }
67 
69 {
70  if (unlikely(std::abs(curvature) < 1.0e-5)) {
71  double sin = (center * axis) / radius;
72  return sin / clamped_sqrt(1 - sqr(sin));
73  } else {
74  double radius2 = sqr(1.0 / curvature);
75  double orthog = clamped_sqrt(radius2 - delta2);
76  Basic2DVector<double> lcenter = center - sgn(curvature) * orthog * axis;
77 
78  double cos = (radius2 + sqr(radius) - lcenter.mag2()) *
79  curvature / (2. * radius);
80  return - cos / clamped_sqrt(1 - sqr(cos));
81  }
82 }
83 
86 {
87  double phi1 = phi(curvature.second, radius);
88  double phi2 = phi(curvature.first, radius);
89 
90  while(unlikely(phi2 < phi1)) phi2 += 2. * M_PI;
91 
92  return Range(phi1 * radius - theTolerance, phi2 * radius + theTolerance);
93 }
94 
96 ThirdHitPredictionFromCircle::curvature(double transverseIP) const
97 {
98  // this is a mess. Use a CAS and lots of drawings to verify...
99 
100  transverseIP = std::abs(transverseIP);
101  double transverseIP2 = sqr(transverseIP);
102  double tip = axis * center;
103  double tip2 = sqr(tip);
104  double lip = axis.x() * center.y() - axis.y() * center.x();
105  double lip2 = sqr(lip);
106 
107  double origin = std::sqrt(tip2 + lip2);
108  double tmp1 = lip2 + tip2 - transverseIP2;
109  double tmp2 = 2. * (tip - transverseIP) * (tip + transverseIP);
110  double tmp3 = 2. * delta * origin;
111  double tmp4 = tmp1 + delta2;
112  double tmp5 = 2. * delta * lip;
113 
114  // I am probably being overly careful here with border cases
115  // but you never know what crap you might get fed
116 
117  double u1, u2;
118  if (unlikely(tmp4 - tmp5 < 1.0e-5)) {
119  u1 = -0.; // yes, I am making use of signed zero
120  u2 = +0.; // -> no -ffast-math please
121  } else {
122  if (unlikely(std::abs(tmp2) < 1.0e-5)) {
123  // the denominator is zero
124  // this means that one of the tracks will be straight
125  // and the other can be computed from the limit of the equation
126  double tmp = lip2 - delta2;
127  u1 = INFINITY; // and require 1 / sqrt(inf^2 + x) = 0 (with x > 0)
128  u2 = (sqr(0.5 * tmp) - delta2 * tip2) / (tmp * tip);
129  if (tip < 0)
130  std::swap(u1, u2);
131  } else {
132  double tmp6 = (tmp4 - tmp5) * (tmp4 + tmp5);
133  if (unlikely(tmp6 < 1.0e-5)) {
134  u1 = -0.;
135  u2 = +0.;
136  } else {
137  double tmp7 = tmp6 > 0 ? (transverseIP * std::sqrt(tmp6) / tmp2) : 0.;
138  double tmp8 = tip * (tmp1 - delta2) / tmp2;
139  // the two quadratic solutions
140  u1 = tmp8 + tmp7;
141  u2 = tmp8 - tmp7;
142  }
143  }
144 
145  if (tmp4 <= std::abs(tmp3)) {
146  if ((tmp3 < 0) == (tip < 0))
147  u2 = +0.;
148  else
149  u1 = -0.;
150  }
151  }
152 
153  return Range(sgn(u1) / std::sqrt(sqr(u1) + delta2),
154  sgn(u2) / std::sqrt(sqr(u2) + delta2));
155 }
156 
158 {
159  Point2D delta = p2 - p1;
160  Point2D axis2 = Point2D(-delta.y(), delta.x()) / delta.mag();
161  Point2D diff = p1 + 0.5 * delta - center;
162  double a = diff.y() * axis2.x() - diff.x() * axis2.y();
163  double b = axis.y() * axis2.x() - axis.x() * axis2.y();
164  return b / a;
165 }
166 
168 {
169  double invDist = invCenterOnAxis(p2);
170  double invDist2 = sqr(invDist);
171  double curv = std::sqrt(invDist2 / (1. + invDist2 * delta2));
172  return sgn(invDist) * curv;
173 }
174 
176 {
177  double invDist = invCenterOnAxis(p2);
178  if (unlikely(std::abs(invDist) < 1.0e-5))
179  return std::abs(p2 * axis);
180  else {
181  double dist = 1.0 / invDist;
182  double radius = std::sqrt(sqr(dist) + delta2);
183  return std::abs((center + axis * dist).mag() - radius);
184  }
185 }
186 
188  const ThirdHitPredictionFromCircle *circle, double z1, double z2, double curv) :
189  circle(circle), curvature(curv), z1(z1)
190 {
191  double absCurv = std::abs(curv);
192  seg = circle->delta;
193 
194  if (likely(absCurv > 1.0e-5)) {
195  seg *= absCurv;
196  seg = seg < -1.0 ? -M_PI_2 : seg > 1.0 ? M_PI_2 : std::asin(seg);
197  seg /= absCurv;
198  }
199 
200  seg *= 2.;
201  dzdu = likely(std::abs(seg) > 1.0e-5) ? ((z2 - z1) / seg) : 99999.0;
202 }
203 
205  const ThirdHitPredictionFromCircle *circle, double z1, double z2, double z3)
206 {
207  static const double maxAngle = M_PI;
208  double halfAngle = (0.5 * maxAngle) * (z2 - z1) / (z3 - z1);
209  if (unlikely(halfAngle <= 0.0))
210  return 0.0;
211 
212  return std::sin(halfAngle) / circle->delta;
213 }
214 
216 {
217  if (unlikely(std::abs(curvature) < 1.0e-5)) {
218  double tip = circle->axis * circle->p1;
219  double lip = circle->axis.y() * circle->p1.x() -
220  circle->axis.x() * circle->p1.y();
221  return z1 + (std::sqrt(sqr(r) - sqr(tip)) - lip) * dzdu;
222  }
223 
224  double radius = 1.0 / curvature;
225  double radius2 = sqr(radius);
226  double orthog = sgn(curvature) * clamped_sqrt(radius2 - circle->delta2);
227  Point2D center = circle->center + orthog * circle->axis;
228 
229  double b2 = center.mag2();
230  double b = std::sqrt(b2);
231 
232  double cos1 = 0.5 * (radius2 + b2 - sqr(r)) * curvature / b;
233  double cos2 = 0.5 * (radius2 + b2 - circle->p1.mag2()) * curvature / b;
234 
235  double phi1 = clamped_acos(cos1);
236  double phi2 = clamped_acos(cos2);
237 
238  // more plausbility checks needed...
239  // the two circles can have two possible intersections
240  double u1 = std::abs((phi1 - phi2) * radius);
241  double u2 = std::abs((phi1 + phi2) * radius);
242 
243  return z1 + ((u1 >= seg && u1 < u2)? u1 : u2) * dzdu;
244 }
245 
247 {
248  if (unlikely(std::abs(dzdu) < 1.0e-5))
249  return 99999.0;
250 
251  if (unlikely(std::abs(curvature) < 1.0e-5)) {
252  double tip = circle->axis * circle->p1;
253  double lip = circle->axis.y() * circle->p1.x() -
254  circle->axis.x() * circle->p1.y();
255  return std::sqrt(sqr(tip) + sqr(lip + (z - z1) / dzdu));
256  }
257 
258  // we won't go below that (see comment below)
259  double minR = (2. * circle->center - circle->p1).mag();
260 
261  double phi = curvature * (z - z1) / dzdu;
262 
263  if (unlikely(std::abs(phi) > 2. * M_PI)) {
264  // with a helix we can get problems here - this is used to get the limits
265  // however, if phi gets large, we get into the regions where we loop back
266  // to smaller r's. The question is - do we care about these tracks?
267  // The answer is probably no: Too much pain, and the rest of the
268  // tracking won't handle looping tracks anyway.
269  // So, what we do here is to return nothing smaller than the radius
270  // than any of the two hits, i.e. the second hit, which is presumably
271  // outside of the 1st hit.
272 
273  return minR;
274  }
275 
276  double radius = 1. / curvature;
277  double orthog = sgn(curvature) * clamped_sqrt(sqr(radius) - circle->delta2);
278  Point2D center = circle->center + orthog * circle->axis;
279  Point2D rel = circle->p1 - center;
280 
281  double c = cos(phi);
282  double s = sin(phi);
283 
284  Point2D p(center.x() + c * rel.x() - s * rel.y(),
285  center.y() + s * rel.x() + c * rel.y());
286 
287  return std::max(minR, p.mag());
288 }
double phi(double curvature, double radius) const
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
float sgn(float val)
Definition: FWPFMaths.cc:9
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T y() const
Definition: PV3DBase.h:62
#define abs(x)
Definition: mlp_lapack.h:159
Range operator()(Range curvature, double radius) const
#define M_PI_2
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
#define unlikely(expr)
double double double z
Range curvature(double transverseIP) const
Basic2DVector< double > Point2D
Basic3DVector< double > Point3D
T curvature(T InversePt, const edm::EventSetup &iSetup)
const T & max(const T &a, const T &b)
double transverseIP(const Basic2DVector< double > &thirdPoint) const
double angle(double curvature, double radius) const
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
T sqrt(T t)
Definition: SSEVec.h:46
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double invCenterOnAxis(const Basic2DVector< double > &thirdPoint) const
double p2[4]
Definition: TauolaWrapper.h:90
T y() const
Cartesian y coordinate.
ThirdHitPredictionFromCircle(const GlobalPoint &P1, const GlobalPoint &P2, float tolerance)
#define M_PI
Definition: BFit3D.cc:3
static double maxCurvature(const ThirdHitPredictionFromCircle *circle, double z1, double z2, double z3)
double b
Definition: hdecay.h:120
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
double p1[4]
Definition: TauolaWrapper.h:89
double a
Definition: hdecay.h:121
Square< F >::type sqr(const F &f)
Definition: Square.h:13
#define likely(expr)
Definition: DDAxes.h:10
Geom::Phi< T > phi() const
long double T
T x() const
Definition: PV3DBase.h:61
T x() const
Cartesian x coordinate.