Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:27:14

0001 #ifndef RecoMuon_TrackerSeedGenerator_L1MuonPixelTrackFitter_H
0002 #define RecoMuon_TrackerSeedGenerator_L1MuonPixelTrackFitter_H
0003 
0004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0005 
0006 #include "DataFormats/GeometryVector/interface/GlobalTag.h"
0007 #include "DataFormats/GeometryVector/interface/Vector3DBase.h"
0008 #include "DataFormats/GeometryVector/interface/Point3DBase.h"
0009 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
0010 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
0011 
0012 #include <vector>
0013 
0014 namespace reco {
0015   class Track;
0016 }
0017 
0018 class TrackingRegion;
0019 class TrackingRecHit;
0020 class L1MuGMTCand;
0021 class PixelRecoLineRZ;
0022 class SeedingHitSet;
0023 class MagneticField;
0024 
0025 class L1MuonPixelTrackFitter {
0026 public:
0027   class Circle {
0028   public:
0029     typedef Vector3DBase<long double, GlobalTag> Vector;
0030     typedef Point3DBase<long double, GlobalTag> Point;
0031     Circle() : theValid(false) {}
0032     Circle(const GlobalPoint& h1, const GlobalPoint& h2, double curvature) : theCurvature(curvature) {
0033       Point p1(h1);
0034       Point p2(h2);
0035       Vector dp = (p2 - p1) / 2.;
0036       int charge = theCurvature > 0 ? 1 : -1;
0037       Vector ec = charge * dp.cross(Vector(0, 0, 1)).unit();
0038       long double dist_tmp = 1. / theCurvature / theCurvature - dp.perp2();
0039       theValid = (dist_tmp > 0.);
0040       theCenter = p1 + dp + ec * sqrt(std::abs(dist_tmp));
0041     }
0042     bool isValid() const { return theValid; }
0043     const Point& center() const { return theCenter; }
0044     const long double& curvature() const { return theCurvature; }
0045 
0046   private:
0047     bool theValid;
0048     long double theCurvature;
0049     Point theCenter;
0050   };
0051 
0052 public:
0053   L1MuonPixelTrackFitter(const edm::ParameterSet& cfg);
0054 
0055   virtual ~L1MuonPixelTrackFitter() {}
0056 
0057   void setL1Constraint(const L1MuGMTCand& muon);
0058   void setPxConstraint(const SeedingHitSet& hits);
0059 
0060   virtual reco::Track* run(const MagneticField& field,
0061                            const std::vector<const TrackingRecHit*>& hits,
0062                            const TrackingRegion& region) const;
0063 
0064   static double getBending(double invPt, double eta, int charge);
0065   static double getBendingError(double invPt, double eta);
0066 
0067 private:
0068   double valInversePt(double phi0, double phiL1, double eta) const;
0069   double errInversePt(double invPt, double eta) const;
0070 
0071   double valPhi(const Circle& c, int charge) const;
0072   double errPhi(double invPt, double eta) const;
0073 
0074   double valCotTheta(const PixelRecoLineRZ& line) const;
0075   double errCotTheta(double invPt, double eta) const;
0076 
0077   double valZip(double curvature, const GlobalPoint& p0, const GlobalPoint& p1) const;
0078   double errZip(double invPt, double eta) const;
0079 
0080   double valTip(const Circle& c, double curvature) const;
0081   double errTip(double invPt, double eta) const;
0082 
0083   double findPt(double phi0, double phiL1, double eta, int charge) const;
0084   double deltaPhi(double phi1, double phi2) const;
0085   static void param(double eta, double& p1, double& p2, double& p3);
0086 
0087 private:
0088   edm::ParameterSet theConfig;
0089 
0090   const double invPtErrorScale;
0091   const double phiErrorScale;
0092   const double cotThetaErrorScale;
0093   const double tipErrorScale;
0094   const double zipErrorScale;
0095 
0096   // L1 constraint
0097   double thePhiL1, theEtaL1;
0098   int theChargeL1;
0099 
0100   // Px constraint
0101   GlobalPoint theHit1, theHit2;
0102 
0103 private:
0104   friend class L1Seeding;
0105 };
0106 #endif