CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PixelFitterByHelixProjections.cc
Go to the documentation of this file.
3 
8 
12 
16 //#include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
18 
23 
26 
30 
31 
33 
34 #include "RZLine.h"
35 #include "CircleFromThreePoints.h"
36 #include "PixelTrackBuilder.h"
38 
39 using namespace std;
40 
41 namespace {
42 
43  int charge(const std::vector<GlobalPoint> & points) {
44  // the cross product will tell me...
45  float dir = (points[1].x()-points[0].x())*(points[2].y()-points[1].y())
46  - (points[1].y()-points[0].y())*(points[2].x()-points[1].x());
47 
48  /*
49  GlobalVector v21 = points[1]-points[0];
50  GlobalVector v32 = points[2]-points[1];
51  float dphi = v32.phi() - v21.phi();
52  while (dphi > Geom::fpi()) dphi -= Geom::ftwoPi();
53  while (dphi < -Geom::fpi()) dphi += Geom::ftwoPi();
54  return (dphi > 0) ? -1 : 1;
55  */
56  return (dir>0) ? -1 : 1;
57  }
58 
59  float cotTheta(const GlobalPoint& inner, const GlobalPoint& outer) {
60  float dr = outer.perp()-inner.perp();
61  float dz = outer.z()-inner.z();
62  return (std::abs(dr) > 1.e-3f) ? dz/dr : 0;
63  }
64 
65  inline float phi(float xC, float yC, int charge) {
66  return (charge>0) ? std::atan2(xC,-yC) : std::atan2(-xC,yC);
67  }
68 
69  float zip(float d0, float phi_p, float curv,
70  const GlobalPoint& pinner, const GlobalPoint& pouter) {
71  //
72  //phi = asin(r*rho/2) with asin(x) ~= x+x**3/(2*3) = x(1+x*x/6);
73  //
74 
75  float phi0 = phi_p - Geom::fhalfPi();
76  GlobalPoint pca(d0*std::cos(phi0), d0*std::sin(phi0),0.);
77 
78  float rho2 = curv*curv;
79  float r1s = (pinner-pca).perp2();
80  double phi1 = std::sqrt(r1s)*(curv*0.5f)*(1.f+r1s*(rho2/24.f));
81  float r2s = (pouter-pca).perp2();
82  double phi2 = std::sqrt(r2s)*(curv*0.5f)*(1.f+r2s*(rho2/24.f));
83  double z1 = pinner.z();
84  double z2 = pouter.z();
85 
86  if (fabs(curv)>1.e-5)
87  return z1 - phi1/(phi1-phi2)*(z1-z2);
88  else {
89  double dr = std::max(std::sqrt(r2s)-std::sqrt(r1s),1.e-5f);
90  return z1-std::sqrt(r1s)*(z2-z1)/dr;
91  }
92  }
93 
94  double errZip( float apt, float eta) {
95  double ziperr=0;
96  double pt = (apt <= 10.) ? apt: 10.;
97  double p1=0, p2=0,p3=0,p4=0;
98  float feta = std::abs(eta);
99  if (feta<=0.8){
100  p1 = 0.12676e-1;
101  p2 = -0.22411e-2;
102  p3 = 0.2987e-3;
103  p4 = -0.12779e-4;
104  } else if (feta <=1.6){
105  p1 = 0.24047e-1;
106  p2 = -0.66935e-2;
107  p3 = 0.88111e-3;
108  p4 = -0.38482e-4;
109  } else {
110  p1 = 0.56084e-1;
111  p2 = -0.13960e-1;
112  p3 = 0.15744e-2;
113  p4 = -0.60757e-4;
114  }
115  ziperr = p1 + p2*pt + p3*pt*pt +p4*pt*pt*pt;
116  return ziperr;
117  }
118 
119  double errZip2( float apt, float eta) {
120  double err = errZip(apt,eta);
121  return err*err;
122  }
123 
124 
125  double errTip(float apt, float eta) {
126  double pt = (apt <= 10.) ? apt : 10.;
127  double p1=0, p2=0;
128  float feta = std::abs(eta);
129  if (feta<=0.8)
130  {
131  p1=5.9e-3;
132  p2=4.7e-3;
133  }
134  else if (feta <=1.6){
135  p1 = 4.9e-3;
136  p2 = 7.1e-3;
137  }
138  else {
139  p1 = 6.4e-3;
140  p2 = 1.0e-2;
141  }
142  double err=0;
143  if (pt != 0) err = (p1 + p2/pt);
144  return err;
145  }
146 
147  double errTip2(float apt, float eta) {
148  double err = errTip(apt,eta);
149  return err*err;
150  }
151 
152 
153 }
154 
155 
157  const edm::ParameterSet& cfg)
158  : theConfig(cfg), theTracker(0), theField(0), theTTRecHitBuilder(0) { }
159 
161  const edm::EventSetup& es,
162  const std::vector<const TrackingRecHit * > & hits,
163  const TrackingRegion & region) const
164 {
165  int nhits = hits.size();
166  if (nhits <2) return 0;
167 
168  vector<GlobalPoint> points(nhits);
169  vector<GlobalError> errors(nhits);
170  vector<bool> isBarrel(nhits);
171 
172  static edm::ESWatcher<TrackerDigiGeometryRecord> watcherTrackerDigiGeometryRecord;
173  if (!theTracker || watcherTrackerDigiGeometryRecord.check(es)) {
175  es.get<TrackerDigiGeometryRecord>().get(trackerESH);
176  theTracker = trackerESH.product();
177  }
178 
179  static edm::ESWatcher<IdealMagneticFieldRecord> watcherIdealMagneticFieldRecord;
180  if (!theField || watcherIdealMagneticFieldRecord.check(es)) {
182  es.get<IdealMagneticFieldRecord>().get(fieldESH);
183  theField = fieldESH.product();
184  }
185 
186  static edm::ESWatcher<TransientRecHitRecord> watcherTransientRecHitRecord;
187  if (!theTTRecHitBuilder || watcherTransientRecHitRecord.check(es)) {
189  std::string builderName = theConfig.getParameter<std::string>("TTRHBuilder");
190  es.get<TransientRecHitRecord>().get(builderName,ttrhbESH);
191  theTTRecHitBuilder = ttrhbESH.product();
192  }
193 
194 
195  for ( int i=0; i!=nhits; ++i) {
197  points[i] = GlobalPoint( recHit->globalPosition().x()-region.origin().x(),
198  recHit->globalPosition().y()-region.origin().y(),
199  recHit->globalPosition().z()-region.origin().z()
200  );
201  errors[i] = recHit->globalPositionError();
202  isBarrel[i] = recHit->detUnit()->type().isBarrel();
203  }
204 
205  CircleFromThreePoints circle = (nhits==2) ?
206  CircleFromThreePoints( GlobalPoint(0.,0.,0.), points[0], points[1]) :
207  CircleFromThreePoints(points[0],points[1],points[2]);
208 
209  float valPhi, valTip, valPt;
210 
211  int iCharge = charge(points);
212  float curvature = circle.curvature();
213 
214  if (curvature > 1.e-4) {
215  float invPt = PixelRecoUtilities::inversePt( circle.curvature(), es);
216  valPt = (invPt > 1.e-4f) ? 1.f/invPt : 1.e4f;
217  CircleFromThreePoints::Vector2D center = circle.center();
218  valTip = iCharge * (center.mag()-1.f/curvature);
219  valPhi = phi(center.x(), center.y(), iCharge);
220  }
221  else {
222  valPt = 1.e4f;
223  GlobalVector direction(points[1]-points[0]);
224  valPhi = direction.phi();
225  valTip = -points[0].x()*sin(valPhi) + points[0].y()*cos(valPhi);
226  }
227 
228  float errPt = 0.055f*valPt + 0.017f*valPt*valPt;
229  float errValTip = errTip(valPt, points.back().eta());
230  float errPhi = 0.002f;
231 
232  float valZip = zip(valTip, valPhi, curvature, points[0],points[1]);
233  float errValZip = errZip(valPt, points.back().eta());
234 
235  float valCotTheta = cotTheta(points[0],points[1]);
236  float errCotTheta = 0.002f;
237 
238  float chi2 = 0;
239  if (nhits > 2) {
240  RZLine rzLine(points,errors,isBarrel);
241  float cottheta, intercept, covss, covii, covsi;
242  rzLine.fit(cottheta, intercept, covss, covii, covsi);
243  chi2 = rzLine.chi2(cottheta, intercept); //FIXME: check which intercept to use!
244  }
245 
246  PixelTrackBuilder builder;
247  Measurement1D pt(valPt, errPt);
248  Measurement1D phi(valPhi, errPhi);
249  Measurement1D cotTheta(valCotTheta, errCotTheta);
250  Measurement1D tip(valTip, errValTip);
251  Measurement1D zip(valZip, errValZip);
252 
253  return builder.build(pt, phi, cotTheta, tip, zip, chi2, iCharge, hits, theField, region.origin() );
254 }
255 
256 
257 
258 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
T perp() const
Definition: PV3DBase.h:66
virtual GlobalPoint origin() const =0
tuple d0
Definition: debug_cff.py:3
PixelFitterByHelixProjections(const edm::ParameterSet &cfg)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Phi< T > phi() const
Definition: PV3DBase.h:63
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
T y() const
Definition: PV3DBase.h:57
#define abs(x)
Definition: mlp_lapack.h:159
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
virtual reco::Track * run(const edm::EventSetup &es, const std::vector< const TrackingRecHit * > &hits, const TrackingRegion &region) const
T eta() const
double charge(const std::vector< uint8_t > &Ampls)
T inversePt(T curvature, const edm::EventSetup &iSetup)
T curvature(T InversePt, const edm::EventSetup &iSetup)
const T & max(const T &a, const T &b)
T sqrt(T t)
Definition: SSEVec.h:28
double p4[4]
Definition: TauolaWrapper.h:92
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
T z() const
Definition: PV3DBase.h:58
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
double f[11][100]
double p2[4]
Definition: TauolaWrapper.h:90
T y() const
Cartesian y coordinate.
Definition: RZLine.h:8
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
T perp2() const
Squared magnitude of transverse component.
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:59
const TransientTrackingRecHitBuilder * theTTRecHitBuilder
double p1[4]
Definition: TauolaWrapper.h:89
float fhalfPi()
Definition: Pi.h:37
void fit(float &cotTheta, float &intercept, float &covss, float &covii, float &covsi) const
Definition: RZLine.cc:44
float chi2(float cotTheta, float intercept) const
Definition: RZLine.cc:50
dbl *** dir
Definition: mlp_gen.cc:35
Definition: DDAxes.h:10
T x() const
Definition: PV3DBase.h:56
T x() const
Cartesian x coordinate.
double p3[4]
Definition: TauolaWrapper.h:91
Definition: DDAxes.h:10