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