CMS 3D CMS Logo

L1MuonPixelTrackFitter.cc
Go to the documentation of this file.
2 
4 
11 
12 #include <cmath>
13 
14 template <class T>
15 T sqr(T t) {
16  return t * t;
17 }
18 
20  : theConfig(cfg),
21  invPtErrorScale{theConfig.getParameter<double>("invPtErrorScale")},
22  phiErrorScale{theConfig.getParameter<double>("phiErrorScale")},
23  cotThetaErrorScale{theConfig.getParameter<double>("cotThetaErrorScale")},
24  tipErrorScale{theConfig.getParameter<double>("tipErrorScale")},
25  zipErrorScale{theConfig.getParameter<double>("zipErrorScale")} {}
26 
28  thePhiL1 = muon.phiValue() + 0.021817;
29  theEtaL1 = muon.etaValue();
30  theChargeL1 = muon.charge();
31 }
32 
34  theHit1 = hits[0]->globalPosition();
35  theHit1 = hits[1]->globalPosition();
36 }
37 
39  const std::vector<const TrackingRecHit*>& hits,
40  const TrackingRegion& region) const {
41  double phi_vtx_fromHits = (theHit2 - theHit1).phi();
42 
43  double invPt = valInversePt(phi_vtx_fromHits, thePhiL1, theEtaL1);
44  double invPtErr = errInversePt(invPt, theEtaL1);
45 
46  int charge = (invPt > 0.) ? 1 : -1;
47  double curvature = PixelRecoUtilities::curvature(invPt, field);
48 
49  Circle circle(theHit1, theHit2, curvature);
50 
51  double valPhi = this->valPhi(circle, charge);
52  double valTip = this->valTip(circle, curvature);
53  double valZip = this->valZip(curvature, theHit1, theHit2);
55 
56  // if ( (fabs(invPt)-0.1)/invPtErr > 3.) return 0;
57 
58  PixelTrackBuilder builder;
59  double valPt = (fabs(invPt) > 1.e-4) ? 1. / fabs(invPt) : 1.e4;
60  double errPt = invPtErrorScale * invPtErr * valPt * valPt;
61  double eta = asinh(valCotTheta);
62  double errPhi = this->errPhi(invPt, eta) * phiErrorScale;
63  double errCotTheta = this->errCotTheta(invPt, eta) * cotThetaErrorScale;
64  double errTip = this->errTip(invPt, eta) * tipErrorScale;
65  double errZip = this->errZip(invPt, eta) * zipErrorScale;
66 
67  Measurement1D pt(valPt, errPt);
72 
73  float chi2 = 0.;
74 
75  /*
76  std::ostringstream str;
77  str <<"\t ipt: " << invPt <<"+/-"<<invPtErr
78  <<"\t pt: " << pt.value() <<"+/-"<<pt.error()
79  <<"\t phi: " << phi.value() <<"+/-"<<phi.error()
80  <<"\t cot: " << cotTheta.value() <<"+/-"<<cotTheta.error()
81  <<"\t tip: " << tip.value() <<"+/-"<<tip.error()
82  <<"\t zip: " << zip.value() <<"+/-"<<zip.error()
83  <<"\t charge: " << charge;
84  std::cout <<str.str()<< std::endl;
85 */
86 
87  return builder.build(pt, phi, cotTheta, tip, zip, chi2, charge, hits, &field);
88 }
89 
90 double L1MuonPixelTrackFitter::valPhi(const Circle& circle, int charge) const {
91  Circle::Point zero(0., 0., 0.);
92  Circle::Vector center = circle.center() - zero;
93  long double radius = center.perp();
94  center /= radius;
95  Circle::Vector dir = center.cross(-charge * Circle::Vector(0., 0., 1.));
96  return dir.phi();
97 }
98 
99 double L1MuonPixelTrackFitter::errPhi(double invPt, double eta) const {
100  //
101  // sigma = p0+p1/|pt|;
102  //
103  double p0, p1;
104  int ieta = int(10 * fabs(eta));
105  switch (ieta) {
106  case 0: {
107  p0 = 0.000597506;
108  p1 = 0.00221057;
109  break;
110  }
111  case 1: {
112  p0 = 0.000591867;
113  p1 = 0.00278744;
114  break;
115  }
116  case 2: {
117  p0 = 0.000635666;
118  p1 = 0.00207433;
119  break;
120  }
121  case 3: {
122  p0 = 0.000619086;
123  p1 = 0.00255121;
124  break;
125  }
126  case 4: {
127  p0 = 0.000572067;
128  p1 = 0.00310618;
129  break;
130  }
131  case 5: {
132  p0 = 0.000596239;
133  p1 = 0.00288442;
134  break;
135  }
136  case 6: {
137  p0 = 0.000607608;
138  p1 = 0.00282996;
139  break;
140  }
141  case 7: {
142  p0 = 0.000606446;
143  p1 = 0.00281118;
144  break;
145  }
146  case 8: {
147  p0 = 0.000594076;
148  p1 = 0.00280546;
149  break;
150  }
151  case 9: {
152  p0 = 0.000579615;
153  p1 = 0.00335534;
154  break;
155  }
156  case 10: {
157  p0 = 0.000659546;
158  p1 = 0.00340443;
159  break;
160  }
161  case 11: {
162  p0 = 0.000659031;
163  p1 = 0.00343151;
164  break;
165  }
166  case 12: {
167  p0 = 0.000738391;
168  p1 = 0.00337297;
169  break;
170  }
171  case 13: {
172  p0 = 0.000798966;
173  p1 = 0.00330008;
174  break;
175  }
176  case 14: {
177  p0 = 0.000702997;
178  p1 = 0.00562643;
179  break;
180  }
181  case 15: {
182  p0 = 0.000973417;
183  p1 = 0.00312666;
184  break;
185  }
186  case 16: {
187  p0 = 0.000995213;
188  p1 = 0.00564278;
189  break;
190  }
191  case 17: {
192  p0 = 0.00121436;
193  p1 = 0.00572704;
194  break;
195  }
196  case 18: {
197  p0 = 0.00119216;
198  p1 = 0.00760204;
199  break;
200  }
201  case 19: {
202  p0 = 0.00141204;
203  p1 = 0.0093777;
204  break;
205  }
206  default: {
207  p0 = 0.00153161;
208  p1 = 0.00940265;
209  break;
210  }
211  }
212  return p0 + p1 * fabs(invPt);
213 }
214 
215 double L1MuonPixelTrackFitter::valCotTheta(const PixelRecoLineRZ& line) const { return line.cotLine(); }
216 
217 double L1MuonPixelTrackFitter::errCotTheta(double invPt, double eta) const {
218  //
219  // sigma = p0+p1/|pt|;
220  //
221  double p0, p1;
222  int ieta = int(10 * fabs(eta));
223  switch (ieta) {
224  case 0: {
225  p0 = 0.00166115;
226  p1 = 5.75533e-05;
227  break;
228  }
229  case 1: {
230  p0 = 0.00157525;
231  p1 = 0.000707437;
232  break;
233  }
234  case 2: {
235  p0 = 0.00122246;
236  p1 = 0.000325456;
237  break;
238  }
239  case 3: {
240  p0 = 0.000852422;
241  p1 = 0.000429216;
242  break;
243  }
244  case 4: {
245  p0 = 0.000637561;
246  p1 = 0.00122298;
247  break;
248  }
249  case 5: {
250  p0 = 0.000555766;
251  p1 = 0.00158096;
252  break;
253  }
254  case 6: {
255  p0 = 0.000641202;
256  p1 = 0.00143339;
257  break;
258  }
259  case 7: {
260  p0 = 0.000803207;
261  p1 = 0.000648816;
262  break;
263  }
264  case 8: {
265  p0 = 0.000741394;
266  p1 = 0.0015289;
267  break;
268  }
269  case 9: {
270  p0 = 0.000652019;
271  p1 = 0.00168873;
272  break;
273  }
274  case 10: {
275  p0 = 0.000716902;
276  p1 = 0.00257556;
277  break;
278  }
279  case 11: {
280  p0 = 0.000800409;
281  p1 = 0.00190563;
282  break;
283  }
284  case 12: {
285  p0 = 0.000808778;
286  p1 = 0.00264139;
287  break;
288  }
289  case 13: {
290  p0 = 0.000775757;
291  p1 = 0.00318478;
292  break;
293  }
294  case 14: {
295  p0 = 0.000705781;
296  p1 = 0.00460576;
297  break;
298  }
299  case 15: {
300  p0 = 0.000580679;
301  p1 = 0.00748248;
302  break;
303  }
304  case 16: {
305  p0 = 0.000561667;
306  p1 = 0.00767487;
307  break;
308  }
309  case 17: {
310  p0 = 0.000521626;
311  p1 = 0.0100178;
312  break;
313  }
314  case 18: {
315  p0 = 0.00064253;
316  p1 = 0.0106062;
317  break;
318  }
319  case 19: {
320  p0 = 0.000636868;
321  p1 = 0.0140047;
322  break;
323  }
324  default: {
325  p0 = 0.000682478;
326  p1 = 0.0163569;
327  break;
328  }
329  }
330  return p0 + p1 * fabs(invPt);
331 }
332 
333 double L1MuonPixelTrackFitter::valTip(const Circle& circle, double curvature) const {
334  Circle::Point zero(0., 0., 0.);
335  Circle::Vector center = circle.center() - zero;
336  long double radius = center.perp();
337  return radius - 1. / fabs(curvature);
338 }
339 
340 double L1MuonPixelTrackFitter::errTip(double invPt, double eta) const {
341  //
342  // sigma = p0+p1/|pt|;
343  //
344  double p0, p1;
345  int ieta = int(10 * fabs(eta));
346  switch (ieta) {
347  case 0: {
348  p0 = 0.00392416;
349  p1 = 0.00551809;
350  break;
351  }
352  case 1: {
353  p0 = 0.00390391;
354  p1 = 0.00543244;
355  break;
356  }
357  case 2: {
358  p0 = 0.0040651;
359  p1 = 0.00406496;
360  break;
361  }
362  case 3: {
363  p0 = 0.00387782;
364  p1 = 0.00797637;
365  break;
366  }
367  case 4: {
368  p0 = 0.00376798;
369  p1 = 0.00866894;
370  break;
371  }
372  case 5: {
373  p0 = 0.0042131;
374  p1 = 0.00462184;
375  break;
376  }
377  case 6: {
378  p0 = 0.00392579;
379  p1 = 0.00784685;
380  break;
381  }
382  case 7: {
383  p0 = 0.00370472;
384  p1 = 0.00790174;
385  break;
386  }
387  case 8: {
388  p0 = 0.00364433;
389  p1 = 0.00928368;
390  break;
391  }
392  case 9: {
393  p0 = 0.00387578;
394  p1 = 0.00640431;
395  break;
396  }
397  case 10: {
398  p0 = 0.00382464;
399  p1 = 0.00960763;
400  break;
401  }
402  case 11: {
403  p0 = 0.0038907;
404  p1 = 0.0104562;
405  break;
406  }
407  case 12: {
408  p0 = 0.00392525;
409  p1 = 0.0106442;
410  break;
411  }
412  case 13: {
413  p0 = 0.00400634;
414  p1 = 0.011218;
415  break;
416  }
417  case 14: {
418  p0 = 0.0036229;
419  p1 = 0.0156403;
420  break;
421  }
422  case 15: {
423  p0 = 0.00444317;
424  p1 = 0.00832987;
425  break;
426  }
427  case 16: {
428  p0 = 0.00465492;
429  p1 = 0.0179908;
430  break;
431  }
432  case 17: {
433  p0 = 0.0049652;
434  p1 = 0.0216647;
435  break;
436  }
437  case 18: {
438  p0 = 0.0051395;
439  p1 = 0.0233692;
440  break;
441  }
442  case 19: {
443  p0 = 0.0062917;
444  p1 = 0.0262175;
445  break;
446  }
447  default: {
448  p0 = 0.00714444;
449  p1 = 0.0253856;
450  break;
451  }
452  }
453  return p0 + p1 * fabs(invPt);
454 }
455 
456 double L1MuonPixelTrackFitter::valZip(double curv, const GlobalPoint& pinner, const GlobalPoint& pouter) const {
457  //
458  // phi = asin(r*rho/2) with asin(x) ~= x+x**3/(2*3)
459  //
460  double rho3 = curv * curv * curv;
461  double r1 = pinner.perp();
462  double phi1 = r1 * curv / 2 + pinner.perp2() * r1 * rho3 / 48.;
463  double r2 = pouter.perp();
464  double phi2 = r2 * curv / 2 + pouter.perp2() * r2 * rho3 / 48.;
465  double z1 = pinner.z();
466  double z2 = pouter.z();
467 
468  return z1 - phi1 / (phi1 - phi2) * (z1 - z2);
469 }
470 
471 double L1MuonPixelTrackFitter::errZip(double invPt, double eta) const {
472  //
473  // sigma = p0+p1/pt;
474  //
475  double p0, p1;
476  int ieta = int(10 * fabs(eta));
477  switch (ieta) {
478  case 0: {
479  p0 = 0.0120743;
480  p1 = 0;
481  break;
482  }
483  case 1: {
484  p0 = 0.0110343;
485  p1 = 0.0051199;
486  break;
487  }
488  case 2: {
489  p0 = 0.00846487;
490  p1 = 0.00570084;
491  break;
492  }
493  case 3: {
494  p0 = 0.00668726;
495  p1 = 0.00331165;
496  break;
497  }
498  case 4: {
499  p0 = 0.00467126;
500  p1 = 0.00578239;
501  break;
502  }
503  case 5: {
504  p0 = 0.0043042;
505  p1 = 0.00598517;
506  break;
507  }
508  case 6: {
509  p0 = 0.00515392;
510  p1 = 0.00495422;
511  break;
512  }
513  case 7: {
514  p0 = 0.0060843;
515  p1 = 0.00320512;
516  break;
517  }
518  case 8: {
519  p0 = 0.00564942;
520  p1 = 0.00478876;
521  break;
522  }
523  case 9: {
524  p0 = 0.00532111;
525  p1 = 0.0073239;
526  break;
527  }
528  case 10: {
529  p0 = 0.00579429;
530  p1 = 0.00952782;
531  break;
532  }
533  case 11: {
534  p0 = 0.00614229;
535  p1 = 0.00977795;
536  break;
537  }
538  case 12: {
539  p0 = 0.00714661;
540  p1 = 0.00550482;
541  break;
542  }
543  case 13: {
544  p0 = 0.0066593;
545  p1 = 0.00999362;
546  break;
547  }
548  case 14: {
549  p0 = 0.00634922;
550  p1 = 0.0148156;
551  break;
552  }
553  case 15: {
554  p0 = 0.00600586;
555  p1 = 0.0318022;
556  break;
557  }
558  case 16: {
559  p0 = 0.00676919;
560  p1 = 0.027456;
561  break;
562  }
563  case 17: {
564  p0 = 0.00670066;
565  p1 = 0.0317005;
566  break;
567  }
568  case 18: {
569  p0 = 0.00752392;
570  p1 = 0.0347714;
571  break;
572  }
573  case 19: {
574  p0 = 0.00791425;
575  p1 = 0.0566665;
576  break;
577  }
578  default: {
579  p0 = 0.00882372;
580  p1 = 0.0596858;
581  break;
582  }
583  }
584  return p0 + p1 * fabs(invPt);
585 }
586 
587 double L1MuonPixelTrackFitter::valInversePt(double phi0, double phiL1, double eta) const {
588  double result = 0.;
589 
590  // solve equtaion p3*result^3 + p1*result + dphi = 0;
591  // where: result = 1/pt
592  // dphi = phiL1 - phi0
593  // Cardan way is used
594 
595  double p1, p2, p3; //parameters p1,p3 are negative by parameter fix in fit!
596  param(eta, p1, p2, p3);
597 
598  double dphi = deltaPhi(phiL1, phi0); // phi0-phiL1
599  if (fabs(dphi) < 0.01) {
600  result = -dphi / p1;
601  } else {
602  double q = dphi / 2. / p3;
603  double p = p1 / 3. / p3; // positive
604  double D = q * q + p * p * p; // positive
605  double u = pow(-q + sqrt(D), 1. / 3.);
606  double v = -pow(q + sqrt(D), 1. / 3.);
607  result = u + v;
608  }
609  return result;
610 }
611 
612 double L1MuonPixelTrackFitter::errInversePt(double invPt, double eta) const {
613  //
614  // pt*sigma(1/pt) = p0+p1*pt;
615  //
616  double p0, p1;
617  int ieta = int(10 * fabs(eta));
618  switch (ieta) {
619  case 0: {
620  p0 = 0.0196835;
621  p1 = 0.00517533;
622  break;
623  }
624  case 1: {
625  p0 = 0.0266583;
626  p1 = 0.00478101;
627  break;
628  }
629  case 2: {
630  p0 = 0.0217164;
631  p1 = 0.00545425;
632  break;
633  }
634  case 3: {
635  p0 = 0.0197547;
636  p1 = 0.00552263;
637  break;
638  }
639  case 4: {
640  p0 = 0.0208778;
641  p1 = 0.00536009;
642  break;
643  }
644  case 5: {
645  p0 = 0.024192;
646  p1 = 0.00521709;
647  break;
648  }
649  case 6: {
650  p0 = 0.0265315;
651  p1 = 0.0051897;
652  break;
653  }
654  case 7: {
655  p0 = 0.0198071;
656  p1 = 0.00566822;
657  break;
658  }
659  case 8: {
660  p0 = 0.0361955;
661  p1 = 0.00486352;
662  break;
663  }
664  case 9: {
665  p0 = 0.037864;
666  p1 = 0.00509094;
667  break;
668  }
669  case 10: {
670  p0 = 0.0382968;
671  p1 = 0.00612354;
672  break;
673  }
674  case 11: {
675  p0 = 0.0308326;
676  p1 = 0.0074234;
677  break;
678  }
679  case 12: {
680  p0 = 0.0248577;
681  p1 = 0.00883049;
682  break;
683  }
684  case 13: {
685  p0 = 0.0279965;
686  p1 = 0.00888293;
687  break;
688  }
689  case 14: {
690  p0 = 0.0372582;
691  p1 = 0.00950252;
692  break;
693  }
694  case 15: {
695  p0 = 0.0281366;
696  p1 = 0.0111501;
697  break;
698  }
699  case 16: {
700  p0 = 0.0421483;
701  p1 = 0.0109413;
702  break;
703  }
704  case 17: {
705  p0 = 0.0461798;
706  p1 = 0.0125824;
707  break;
708  }
709  case 18: {
710  p0 = 0.0530603;
711  p1 = 0.0132638;
712  break;
713  }
714  case 19: {
715  p0 = 0.0601148;
716  p1 = 0.0147911;
717  break;
718  }
719  default: {
720  p0 = 0.0552377;
721  p1 = 0.0155574;
722  break;
723  }
724  }
725  return p1 + p0 * fabs(invPt);
726 }
727 
728 double L1MuonPixelTrackFitter::findPt(double phi0, double phiL1, double eta, int charge) const {
729  double dphi_min = fabs(deltaPhi(phi0, phiL1));
730  double pt_best = 1.;
731  double pt_cur = 1;
732  while (pt_cur < 10000.) {
733  double phi_exp = phi0 + getBending(1. / pt_cur, eta, charge);
734  double dphi = fabs(deltaPhi(phi_exp, phiL1));
735  if (dphi < dphi_min) {
736  pt_best = pt_cur;
737  dphi_min = dphi;
738  }
739  if (pt_cur < 10.)
740  pt_cur += 0.01;
741  else if (pt_cur < 20.)
742  pt_cur += 0.025;
743  else if (pt_cur < 100.)
744  pt_cur += 0.1;
745  else
746  pt_cur += 1;
747  };
748  return pt_best;
749 }
750 
751 double L1MuonPixelTrackFitter::getBending(double invPt, double eta, int charge) {
752  double p1, p2, p3;
753  param(eta, p1, p2, p3);
754  return charge * p1 * invPt + charge * p2 * invPt * invPt + charge * p3 * invPt * invPt * invPt;
755 }
756 
757 double L1MuonPixelTrackFitter::getBendingError(double invPt, double eta) {
758  int ieta = int(10 * fabs(eta));
759  double p0, p1;
760  switch (ieta) {
761  case 0: {
762  p0 = 0.0196835;
763  p1 = 0.00517533;
764  break;
765  }
766  case 1: {
767  p0 = 0.0266583;
768  p1 = 0.00478101;
769  break;
770  }
771  case 2: {
772  p0 = 0.0217164;
773  p1 = 0.00545425;
774  break;
775  }
776  case 3: {
777  p0 = 0.0197547;
778  p1 = 0.00552263;
779  break;
780  }
781  case 4: {
782  p0 = 0.0208778;
783  p1 = 0.00536009;
784  break;
785  }
786  case 5: {
787  p0 = 0.024192;
788  p1 = 0.00521709;
789  break;
790  }
791  case 6: {
792  p0 = 0.0265315;
793  p1 = 0.0051897;
794  break;
795  }
796  case 7: {
797  p0 = 0.0198071;
798  p1 = 0.00566822;
799  break;
800  }
801  case 8: {
802  p0 = 0.0361955;
803  p1 = 0.00486352;
804  break;
805  }
806  case 9: {
807  p0 = 0.037864;
808  p1 = 0.00509094;
809  break;
810  }
811  case 10: {
812  p0 = 0.0382968;
813  p1 = 0.00612354;
814  break;
815  }
816  case 11: {
817  p0 = 0.0308326;
818  p1 = 0.0074234;
819  break;
820  }
821  case 12: {
822  p0 = 0.0248577;
823  p1 = 0.00883049;
824  break;
825  }
826  case 13: {
827  p0 = 0.0279965;
828  p1 = 0.00888293;
829  break;
830  }
831  case 14: {
832  p0 = 0.0372582;
833  p1 = 0.00950252;
834  break;
835  }
836  case 15: {
837  p0 = 0.0281366;
838  p1 = 0.0111501;
839  break;
840  }
841  case 16: {
842  p0 = 0.0421483;
843  p1 = 0.0109413;
844  break;
845  }
846  case 17: {
847  p0 = 0.0461798;
848  p1 = 0.0125824;
849  break;
850  }
851  case 18: {
852  p0 = 0.0530603;
853  p1 = 0.0132638;
854  break;
855  }
856  case 19: {
857  p0 = 0.0601148;
858  p1 = 0.0147911;
859  break;
860  }
861  default: {
862  p0 = 0.0552377;
863  p1 = 0.0155574;
864  break;
865  }
866  }
867  return p0 + p1 * sqr(invPt);
868 }
869 
870 void L1MuonPixelTrackFitter::param(double eta, double& p1, double& p2, double& p3) {
871  int ieta = int(10 * fabs(eta));
872  switch (ieta) {
873  case 0: {
874  p1 = -2.68016;
875  p2 = 0;
876  p3 = -12.9653;
877  break;
878  }
879  case 1: {
880  p1 = -2.67864;
881  p2 = 0;
882  p3 = -12.0036;
883  break;
884  }
885  case 2: {
886  p1 = -2.72997;
887  p2 = 0;
888  p3 = -10.3468;
889  break;
890  }
891  case 3: {
892  p1 = -2.68836;
893  p2 = 0;
894  p3 = -12.3369;
895  break;
896  }
897  case 4: {
898  p1 = -2.66885;
899  p2 = 0;
900  p3 = -11.589;
901  break;
902  }
903  case 5: {
904  p1 = -2.64932;
905  p2 = 0;
906  p3 = -12.7176;
907  break;
908  }
909  case 6: {
910  p1 = -2.64513;
911  p2 = 0;
912  p3 = -11.3912;
913  break;
914  }
915  case 7: {
916  p1 = -2.64487;
917  p2 = 0;
918  p3 = -9.83239;
919  break;
920  }
921  case 8: {
922  p1 = -2.58875;
923  p2 = 0;
924  p3 = -10.0541;
925  break;
926  }
927  case 9: {
928  p1 = -2.5551;
929  p2 = 0;
930  p3 = -9.75757;
931  break;
932  }
933  case 10: {
934  p1 = -2.55294;
935  p2 = 0;
936  p3 = -7.25238;
937  break;
938  }
939  case 11: {
940  p1 = -2.45333;
941  p2 = 0;
942  p3 = -7.966;
943  break;
944  }
945  case 12: {
946  p1 = -2.29283;
947  p2 = 0;
948  p3 = -7.03231;
949  break;
950  }
951  case 13: {
952  p1 = -2.13167;
953  p2 = 0;
954  p3 = -4.29182;
955  break;
956  }
957  case 14: {
958  p1 = -1.9102;
959  p2 = 0;
960  p3 = -3.21295;
961  break;
962  }
963  case 15: {
964  p1 = -1.74552;
965  p2 = 0;
966  p3 = -0.531827;
967  break;
968  }
969  case 16: {
970  p1 = -1.53642;
971  p2 = 0;
972  p3 = -1.74057;
973  break;
974  }
975  case 17: {
976  p1 = -1.39446;
977  p2 = 0;
978  p3 = -1.56819;
979  break;
980  }
981  case 18: {
982  p1 = -1.26176;
983  p2 = 0;
984  p3 = -0.598631;
985  break;
986  }
987  case 19: {
988  p1 = -1.14133;
989  p2 = 0;
990  p3 = -0.0941055;
991  break;
992  }
993  default: {
994  p1 = -1.02086;
995  p2 = 0;
996  p3 = -0.100491;
997  break;
998  }
999  }
1000 }
1001 
1002 double L1MuonPixelTrackFitter::deltaPhi(double phi1, double phi2) const {
1003  while (phi1 >= 2 * M_PI)
1004  phi1 -= 2 * M_PI;
1005  while (phi2 >= 2 * M_PI)
1006  phi2 -= 2 * M_PI;
1007  while (phi1 < 0)
1008  phi1 += 2 * M_PI;
1009  while (phi2 < 0)
1010  phi2 += 2 * M_PI;
1011  double dPhi = phi2 - phi1;
1012 
1013  if (dPhi > M_PI)
1014  dPhi -= 2 * M_PI;
1015  if (dPhi < -M_PI)
1016  dPhi += 2 * M_PI;
1017 
1018  return dPhi;
1019 }
double valInversePt(double phi0, double phiL1, double eta) const
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
T perp() const
Definition: PV3DBase.h:69
Vector3DBase< typename PreciseFloatType< T, U >::Type, FrameTag > cross(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:110
void setPxConstraint(const SeedingHitSet &hits)
T z() const
Definition: PV3DBase.h:61
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float zip(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:90
double errInversePt(double invPt, double eta) const
static double getBendingError(double invPt, double eta)
T curvature(T InversePt, const MagneticField &field)
reco::Track * build(const Measurement1D &pt, const Measurement1D &phi, const Measurement1D &cotTheta, const Measurement1D &tip, const Measurement1D &zip, float chi2, int charge, const std::vector< const TrackingRecHit *> &hits, const MagneticField *mf, const GlobalPoint &reference=GlobalPoint(0, 0, 0)) const
T sqr(T t)
double errTip(double invPt, double eta) const
static void param(double eta, double &p1, double &p2, double &p3)
T sqrt(T t)
Definition: SSEVec.h:23
double valTip(const Circle &c, double curvature) const
virtual reco::Track * run(const MagneticField &field, const std::vector< const TrackingRecHit *> &hits, const TrackingRegion &region) const
static double getBending(double invPt, double eta, int charge)
double errCotTheta(double invPt, double eta) const
double deltaPhi(double phi1, double phi2) const
#define M_PI
double valZip(double curvature, const GlobalPoint &p0, const GlobalPoint &p1) const
DecomposeProduct< arg, typename Div::arg > D
Definition: Factorize.h:141
T perp2() const
Definition: PV3DBase.h:68
double errPhi(double invPt, double eta) const
double errZip(double invPt, double eta) const
double findPt(double phi0, double phiL1, double eta, int charge) const
void setL1Constraint(const L1MuGMTCand &muon)
long double T
L1MuonPixelTrackFitter(const edm::ParameterSet &cfg)
double valPhi(const Circle &c, int charge) const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
double valCotTheta(const PixelRecoLineRZ &line) const