CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackFitter.cc
Go to the documentation of this file.
3 
8 
12 
17 
22 
25 
28 
30 
34 
35 using namespace std;
36 
37 /*****************************************************************************/
39  (const edm::ParameterSet& cfg) :
40  theConfig(cfg), theTracker(0), theField(0), theTTRecHitBuilder(0)
41 {
42 }
43 
44 /*****************************************************************************/
46  (const edm::EventSetup& es,
47  const std::vector<const TrackingRecHit *> & hits,
48  const TrackingRegion & region) const
49 {
50  int nhits = hits.size();
51  if(nhits <2) return 0;
52 
53  vector<GlobalPoint> points;
54  vector<GlobalError> errors;
55  vector<bool> isBarrel;
56 
57  if (!theField || !theTracker || !theTTRecHitBuilder)
58  {
60  es.get<TrackerDigiGeometryRecord>().get(trackerESH);
61  theTracker = trackerESH.product();
62 
64  es.get<IdealMagneticFieldRecord>().get(fieldESH);
65  theField = fieldESH.product();
66 
68  std::string builderName = theConfig.getParameter<std::string>("TTRHBuilder");
69  es.get<TransientRecHitRecord>().get(builderName,ttrhbESH);
70  theTTRecHitBuilder = ttrhbESH.product();
71  }
72 
73  for (vector<const TrackingRecHit*>::const_iterator ih = hits.begin();
74  ih!= hits.end(); ih++)
75  {
77  theTTRecHitBuilder->build(*ih);
78 
79  points.push_back(recHit->globalPosition());
80  errors.push_back(recHit->globalPositionError());
81  isBarrel.push_back(recHit->detUnit()->type().isBarrel() );
82  }
83 
84  CircleFromThreePoints circle = (nhits==2) ?
85  CircleFromThreePoints(GlobalPoint(0.,0.,0.), points[0], points[1]) :
86  CircleFromThreePoints(points[0],points[1],points[2]);
87 
88  int charge = getCharge(points);
89  float curvature = circle.curvature();
90 
91  // pt
92  float invPt = PixelRecoUtilities::inversePt(curvature, es);
93  float valPt = (invPt > 1.e-4) ? 1./invPt : 1.e4;
94  float errPt = 0.055*valPt + 0.017*valPt*valPt;
95 
96  CircleFromThreePoints::Vector2D center = circle.center();
97 
98  // tip
99  float valTip = charge * (center.mag()-1/curvature);
100  // zip
101  float valZip = getZip(valTip, curvature, points[0],points[1]);
102  // phi
103  float valPhi = getPhi(center.x(), center.y(), charge);
104  // cot(theta), update zip
105  float valCotTheta =
106  getCotThetaAndUpdateZip(points[0],points[1], 1/curvature,
107  valPhi,valTip,valZip);
108 
109  // errors
110  float errTip, errZip;
111  getErrTipAndErrZip(valPt, points.back().eta(), errTip,errZip);
112  float errPhi = 0.002;
113  float errCotTheta = 0.002;
114 
115  float chi2 = 0;
116  if(nhits > 2)
117  {
118  RZLine rzLine(points,errors,isBarrel);
119  float cotTheta, intercept, covss, covii, covsi;
120  rzLine.fit(cotTheta, intercept, covss, covii, covsi);
121  chi2 = rzLine.chi2(cotTheta, intercept);
122  //FIXME: check which intercept to use!
123  }
124 
125  // build pixel track
126  PixelTrackBuilder builder;
127 
128  Measurement1D pt (valPt, errPt);
129  Measurement1D phi (valPhi, errPhi);
130  Measurement1D cotTheta(valCotTheta, errCotTheta);
131  Measurement1D tip (valTip, errTip);
132  Measurement1D zip (valZip, errZip);
133 
134  return builder.build(pt, phi, cotTheta, tip, zip, chi2,
135  charge, hits, theField);
136 }
137 
138 /*****************************************************************************/
140  (const vector<GlobalPoint> & points) const
141 {
142  GlobalVector v21 = points[1]-points[0];
143  GlobalVector v32 = points[2]-points[1];
144  float dphi = v32.phi() - v21.phi();
145  while (dphi > M_PI) dphi -= 2*M_PI;
146  while (dphi < -M_PI) dphi += 2*M_PI;
147  return (dphi > 0) ? -1 : 1;
148 }
149 
150 /*****************************************************************************/
153  float radius, float phi, float d0, float& zip) const
154 {
155  float chi = phi - M_PI_2;
156  GlobalPoint IP(d0*cos(chi), d0*sin(chi),zip);
157 
158  float phi1 = 2*asin(0.5*(inner - IP).perp()/radius);
159  float phi2 = 2*asin(0.5*(outer - IP).perp()/radius);
160 
161  float dr = radius*(phi2 - phi1);
162  float dz = outer.z()-inner.z();
163 
164  // Recalculate ZIP
165  zip = (inner.z()*phi2 - outer.z()*phi1)/(phi2 - phi1);
166 
167  return (fabs(dr) > 1.e-3) ? dz/dr : 0;
168 }
169 
170 /*****************************************************************************/
172  (float xC, float yC, int charge) const
173 {
174  float phiC;
175 
176  if (charge>0) phiC = atan2(xC,-yC);
177  else phiC = atan2(-xC,yC);
178 
179  return phiC;
180 }
181 
182 /*****************************************************************************/
184  (float d0, float curv,
185  const GlobalPoint& inner, const GlobalPoint& outer) const
186 {
187  // phi = asin(r*rho/2) with asin(x) ~= x+x**3/(2*3)
188  float rho3 = curv*curv*curv;
189 
190  float r1 = inner.perp();
191  double phi1 = r1*curv/2 + inner.perp2()*r1*rho3/48.;
192 
193  float r2 = outer.perp();
194  double phi2 = r2*curv/2 + outer.perp2()*r2*rho3/48.;
195 
196  double z1 = inner.z();
197  double z2 = outer.z();
198 
199  return z1 - phi1/(phi1-phi2)*(z1-z2);
200 }
201 
202 /*****************************************************************************/
204  (float pt, float eta, float & errTip, float & errZip) const
205 {
206  float coshEta = cosh(eta);
207 
208  { // transverse
209  float c_ms = 0.0115; //0.0115;
210  float s_le = 0.0095; //0.0123;
211  float s_ms2 = c_ms*c_ms / (pt*pt) * coshEta;
212 
213  errTip = sqrt(s_le*s_le + s_ms2 );
214  }
215 
216  { // z
217  float c_ms = 0.0070;
218  float s_le = 0.0135;
219 
220  errZip = sqrt( (s_le*s_le + c_ms*c_ms/(pt*pt)) * coshEta*coshEta*coshEta);
221  }
222 }
223 
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:204
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
T perp2() const
Definition: PV3DBase.h:71
T eta() const
int getCharge(const std::vector< GlobalPoint > &points) const
Definition: TrackFitter.cc:140
double charge(const std::vector< uint8_t > &Ampls)
T inversePt(T curvature, const edm::EventSetup &iSetup)
float getZip(float d0, float curv, const GlobalPoint &inner, const GlobalPoint &outer) const
Definition: TrackFitter.cc:184
T curvature(T InversePt, const edm::EventSetup &iSetup)
T sqrt(T t)
Definition: SSEVec.h:48
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
tuple IP
Definition: listHistos.py:76
float getCotThetaAndUpdateZip(const GlobalPoint &inner, const GlobalPoint &outer, float radius, float phi, float d0, float &zip) const
Definition: TrackFitter.cc:152
T y() const
Cartesian y coordinate.
Definition: RZLine.h:8
#define M_PI
Definition: BFit3D.cc:3
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
T perp() const
Magnitude of transverse component.
void fit(float &cotTheta, float &intercept, float &covss, float &covii, float &covsi) const
Definition: RZLine.cc:44
TrackFitter(const edm::ParameterSet &cfg)
Definition: TrackFitter.cc:39
float chi2(float cotTheta, float intercept) const
Definition: RZLine.cc:50
float getPhi(float xC, float yC, int charge) const
Definition: TrackFitter.cc:172
T x() const
Cartesian x coordinate.
virtual reco::Track * run(const edm::EventSetup &es, const std::vector< const TrackingRecHit * > &hits, const TrackingRegion &region) const
Definition: TrackFitter.cc:46
Definition: DDAxes.h:10