CMS 3D CMS Logo

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