CMS 3D CMS Logo

L1MuonPixelTrackFitter.h
Go to the documentation of this file.
1 #ifndef RecoMuon_TrackerSeedGenerator_L1MuonPixelTrackFitter_H
2 #define RecoMuon_TrackerSeedGenerator_L1MuonPixelTrackFitter_H
3 
5 
11 
12 #include <vector>
13 
14 namespace edm {
15  class EventSetup;
16 }
17 namespace reco {
18  class Track;
19 }
20 
21 class TrackingRegion;
22 class TrackingRecHit;
23 class L1MuGMTCand;
24 class PixelRecoLineRZ;
25 class SeedingHitSet;
26 
28 public:
29  class Circle {
30  public:
34  Circle(const GlobalPoint& h1, const GlobalPoint& h2, double curvature) : theCurvature(curvature) {
35  Point p1(h1);
36  Point p2(h2);
37  Vector dp = (p2 - p1) / 2.;
38  int charge = theCurvature > 0 ? 1 : -1;
39  Vector ec = charge * dp.cross(Vector(0, 0, 1)).unit();
40  long double dist_tmp = 1. / theCurvature / theCurvature - dp.perp2();
41  theValid = (dist_tmp > 0.);
42  theCenter = p1 + dp + ec * sqrt(std::abs(dist_tmp));
43  }
44  bool isValid() const { return theValid; }
45  const Point& center() const { return theCenter; }
46  const long double& curvature() const { return theCurvature; }
47 
48  private:
49  bool theValid;
50  long double theCurvature;
52  };
53 
54 public:
56 
58 
59  void setL1Constraint(const L1MuGMTCand& muon);
60  void setPxConstraint(const SeedingHitSet& hits);
61 
62  virtual reco::Track* run(const edm::EventSetup& es,
63  const std::vector<const TrackingRecHit*>& hits,
64  const TrackingRegion& region) const;
65 
66  static double getBending(double invPt, double eta, int charge);
67  static double getBendingError(double invPt, double eta);
68 
69 private:
70  double valInversePt(double phi0, double phiL1, double eta) const;
71  double errInversePt(double invPt, double eta) const;
72 
73  double valPhi(const Circle& c, int charge) const;
74  double errPhi(double invPt, double eta) const;
75 
76  double valCotTheta(const PixelRecoLineRZ& line) const;
77  double errCotTheta(double invPt, double eta) const;
78 
79  double valZip(double curvature, const GlobalPoint& p0, const GlobalPoint& p1) const;
80  double errZip(double invPt, double eta) const;
81 
82  double valTip(const Circle& c, double curvature) const;
83  double errTip(double invPt, double eta) const;
84 
85  double findPt(double phi0, double phiL1, double eta, int charge) const;
86  double deltaPhi(double phi1, double phi2) const;
87  static void param(double eta, double& p1, double& p2, double& p3);
88 
89 private:
91 
92  const double invPtErrorScale;
93  const double phiErrorScale;
94  const double cotThetaErrorScale;
95  const double tipErrorScale;
96  const double zipErrorScale;
97 
98  // L1 constraint
99  double thePhiL1, theEtaL1;
101 
102  // Px constraint
104 
105 private:
106  friend class L1Seeding;
107 };
108 #endif
double valInversePt(double phi0, double phiL1, double eta) const
Vector3DBase< long double, GlobalTag > Vector
Point3DBase< long double, GlobalTag > Point
void setPxConstraint(const SeedingHitSet &hits)
const long double & curvature() const
double errInversePt(double invPt, double eta) const
static double getBendingError(double invPt, double eta)
T curvature(T InversePt, const MagneticField &field)
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:19
double valTip(const Circle &c, double curvature) const
virtual reco::Track * run(const edm::EventSetup &es, const std::vector< const TrackingRecHit *> &hits, const TrackingRegion &region) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static double getBending(double invPt, double eta, int charge)
double errCotTheta(double invPt, double eta) const
double deltaPhi(double phi1, double phi2) const
double valZip(double curvature, const GlobalPoint &p0, const GlobalPoint &p1) const
fixed size matrix
HLT enums.
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)
Circle(const GlobalPoint &h1, const GlobalPoint &h2, double curvature)
L1MuonPixelTrackFitter(const edm::ParameterSet &cfg)
double valPhi(const Circle &c, int charge) const
double valCotTheta(const PixelRecoLineRZ &line) const