CMS 3D CMS Logo

TrackFitter.cc
Go to the documentation of this file.
22 
23 using namespace std;
24 
25 namespace {
26 
28  GlobalVector v21 = points[1] - points[0];
29  GlobalVector v32 = points[2] - points[1];
30  float dphi = v32.phi() - v21.phi();
31  while (dphi > M_PI)
32  dphi -= 2 * M_PI;
33  while (dphi < -M_PI)
34  dphi += 2 * M_PI;
35  return (dphi > 0) ? -1 : 1;
36  }
37 
38 } // namespace
39 
40 /*****************************************************************************/
41 std::unique_ptr<reco::Track> TrackFitter::run(const std::vector<const TrackingRecHit*>& hits,
42  const TrackingRegion& region,
43  const edm::EventSetup& setup) const {
44  std::unique_ptr<reco::Track> ret;
45 
46  int nhits = hits.size();
47  if (nhits < 2)
48  return ret;
49 
50  declareDynArray(GlobalPoint, nhits, points);
52  declareDynArray(bool, nhits, isBarrel);
53 
54  unsigned int i = 0;
55  for (auto const& ih : hits) {
56  auto recHit = theTTRecHitBuilder->build(ih);
57 
58  points[i] = recHit->globalPosition();
59  errors[i] = recHit->globalPositionError();
60  isBarrel[++i] = recHit->detUnit()->type().isBarrel();
61  }
62 
63  CircleFromThreePoints circle = (nhits == 2) ? CircleFromThreePoints(GlobalPoint(0., 0., 0.), points[0], points[1])
64  : CircleFromThreePoints(points[0], points[1], points[2]);
65 
66  int charge = getCharge(points);
67  float curvature = circle.curvature();
68 
69  // pt
70  float invPt = PixelRecoUtilities::inversePt(curvature, setup);
71  float valPt = (invPt > 1.e-4) ? 1. / invPt : 1.e4;
72  float errPt = 0.055 * valPt + 0.017 * valPt * valPt;
73 
74  CircleFromThreePoints::Vector2D center = circle.center();
75 
76  // tip
77  float valTip = charge * (center.mag() - 1 / curvature);
78  // zip
79  float valZip = getZip(valTip, curvature, points[0], points[1]);
80  // phi
81  float valPhi = getPhi(center.x(), center.y(), charge);
82  // cot(theta), update zip
83  float valCotTheta = getCotThetaAndUpdateZip(points[0], points[1], 1 / curvature, valPhi, valTip, valZip);
84 
85  // errors
86  float errTip, errZip;
87  getErrTipAndErrZip(valPt, points.back().eta(), errTip, errZip);
88  float errPhi = 0.002;
89  float errCotTheta = 0.002;
90 
91  float chi2 = 0;
92  if (nhits > 2) {
93  RZLine rzLine(points, errors, isBarrel);
94  chi2 = rzLine.chi2();
95  }
96 
97  // build pixel track
98  PixelTrackBuilder builder;
99 
100  Measurement1D pt(valPt, errPt);
101  Measurement1D phi(valPhi, errPhi);
102  Measurement1D cotTheta(valCotTheta, errCotTheta);
103  Measurement1D tip(valTip, errTip);
104  Measurement1D zip(valZip, errZip);
105 
106  ret.reset(builder.build(pt, phi, cotTheta, tip, zip, chi2, charge, hits, theField));
107  return ret;
108 }
109 
110 /*****************************************************************************/
112  const GlobalPoint& inner, const GlobalPoint& outer, float radius, float phi, float d0, float& zip) const {
113  float chi = phi - M_PI_2;
114  GlobalPoint IP(d0 * cos(chi), d0 * sin(chi), zip);
115 
116  float phi1 = 2 * asin(0.5 * (inner - IP).perp() / radius);
117  float phi2 = 2 * asin(0.5 * (outer - IP).perp() / radius);
118 
119  float dr = radius * (phi2 - phi1);
120  float dz = outer.z() - inner.z();
121 
122  // Recalculate ZIP
123  zip = (inner.z() * phi2 - outer.z() * phi1) / (phi2 - phi1);
124 
125  return (fabs(dr) > 1.e-3) ? dz / dr : 0;
126 }
127 
128 /*****************************************************************************/
129 float TrackFitter::getPhi(float xC, float yC, int charge) const {
130  float phiC;
131 
132  if (charge > 0)
133  phiC = atan2(xC, -yC);
134  else
135  phiC = atan2(-xC, yC);
136 
137  return phiC;
138 }
139 
140 /*****************************************************************************/
141 float TrackFitter::getZip(float d0, float curv, const GlobalPoint& inner, const GlobalPoint& outer) const {
142  // phi = asin(r*rho/2) with asin(x) ~= x+x**3/(2*3)
143  float rho3 = curv * curv * curv;
144 
145  float r1 = inner.perp();
146  double phi1 = r1 * curv / 2 + inner.perp2() * r1 * rho3 / 48.;
147 
148  float r2 = outer.perp();
149  double phi2 = r2 * curv / 2 + outer.perp2() * r2 * rho3 / 48.;
150 
151  double z1 = inner.z();
152  double z2 = outer.z();
153 
154  return z1 - phi1 / (phi1 - phi2) * (z1 - z2);
155 }
156 
157 /*****************************************************************************/
158 void TrackFitter::getErrTipAndErrZip(float pt, float eta, float& errTip, float& errZip) const {
159  float coshEta = cosh(eta);
160 
161  { // transverse
162  float c_ms = 0.0115; //0.0115;
163  float s_le = 0.0095; //0.0123;
164  float s_ms2 = c_ms * c_ms / (pt * pt) * coshEta;
165 
166  errTip = sqrt(s_le * s_le + s_ms2);
167  }
168 
169  { // z
170  float c_ms = 0.0070;
171  float s_le = 0.0135;
172 
173  errZip = sqrt((s_le * s_le + c_ms * c_ms / (pt * pt)) * coshEta * coshEta * coshEta);
174  }
175 }
#define declareDynArray(T, n, x)
Definition: DynArray.h:91
int getCharge(const SiStripCluster *cluster, int &nSatStrip, const GeomDetUnit &detUnit, const std::vector< std::vector< float > > &calibGains, const unsigned int &m_off)
Definition: DeDxTools.cc:214
T perp() const
Definition: PV3DBase.h:69
ret
prodAgent to be discontinued
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
#define M_PI_2
void getErrTipAndErrZip(float pt, float eta, float &errZip, float &errTip) const
Definition: TrackFitter.cc:158
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
std::unique_ptr< reco::Track > run(const std::vector< const TrackingRecHit * > &hits, const TrackingRegion &region, const edm::EventSetup &setup) const override
Definition: TrackFitter.cc:41
T perp2() const
Definition: PV3DBase.h:68
T inversePt(T curvature, const edm::EventSetup &iSetup)
float getZip(float d0, float curv, const GlobalPoint &inner, const GlobalPoint &outer) const
Definition: TrackFitter.cc:141
T curvature(T InversePt, const edm::EventSetup &iSetup)
OutputIterator zip(InputIterator1 first1, InputIterator1 last1, InputIterator2 first2, InputIterator2 last2, OutputIterator result, Compare comp)
T sqrt(T t)
Definition: SSEVec.h:19
T z() const
Definition: PV3DBase.h:61
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
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
float getCotThetaAndUpdateZip(const GlobalPoint &inner, const GlobalPoint &outer, float radius, float phi, float d0, float &zip) const
Definition: TrackFitter.cc:111
T y() const
Cartesian y coordinate.
#define M_PI
Definition: RZLine.h:12
T perp() const
Magnitude of transverse component.
Definition: errors.py:1
float getPhi(float xC, float yC, int charge) const
Definition: TrackFitter.cc:129
float chi2() const
Definition: RZLine.h:94
T x() const
Cartesian x coordinate.