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 
29 
30 
32 
33 #include "RZLine.h"
34 #include "CircleFromThreePoints.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 func_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 
96  const edm::ParameterSet& cfg)
97  : theConfig(cfg), theTracker(0), theField(0), theTTRecHitBuilder(0) {}
98 
100  const edm::EventSetup& es,
101  const std::vector<const TrackingRecHit * > & hits,
102  const TrackingRegion & region) const
103 {
104  int nhits = hits.size();
105  if (nhits <2) return 0;
106 
107  vector<GlobalPoint> points(nhits);
108  vector<GlobalError> errors(nhits);
109  vector<bool> isBarrel(nhits);
110 
111  if (theTrackerWatcher.check(es)) {
113  es.get<TrackerDigiGeometryRecord>().get(trackerESH);
114  theTracker = trackerESH.product();
115  }
116 
117  if (theFieldWatcher.check(es)) {
119  es.get<IdealMagneticFieldRecord>().get(fieldESH);
120  theField = fieldESH.product();
121  }
122 
125  std::string builderName = theConfig.getParameter<std::string>("TTRHBuilder");
126  es.get<TransientRecHitRecord>().get(builderName,ttrhbESH);
127  theTTRecHitBuilder = ttrhbESH.product();
128  }
129 
130 
131  for ( int i=0; i!=nhits; ++i) {
132  auto const & recHit = hits[i];
133  points[i] = GlobalPoint( recHit->globalPosition().x()-region.origin().x(),
134  recHit->globalPosition().y()-region.origin().y(),
135  recHit->globalPosition().z()-region.origin().z()
136  );
137  errors[i] = recHit->globalPositionError();
138  isBarrel[i] = recHit->detUnit()->type().isBarrel();
139  }
140 
141  CircleFromThreePoints circle = (nhits==2) ?
142  CircleFromThreePoints( GlobalPoint(0.,0.,0.), points[0], points[1]) :
143  CircleFromThreePoints(points[0],points[1],points[2]);
144 
145  float valPhi, valTip, valPt;
146 
147  int iCharge = charge(points);
148  float curvature = circle.curvature();
149 
150  if ((curvature > 1.e-4)&&
151  (likely(theField->inTesla(GlobalPoint(0.,0.,0.)).z()>0.01))) {
152  float invPt = PixelRecoUtilities::inversePt( circle.curvature(), es);
153  valPt = (invPt > 1.e-4f) ? 1.f/invPt : 1.e4f;
154  CircleFromThreePoints::Vector2D center = circle.center();
155  valTip = iCharge * (center.mag()-1.f/curvature);
156  valPhi = func_phi(center.x(), center.y(), iCharge);
157  }
158  else {
159  valPt = 1.e4f;
160  GlobalVector direction(points[1]-points[0]);
161  valPhi = direction.phi();
162  valTip = -points[0].x()*sin(valPhi) + points[0].y()*cos(valPhi);
163  }
164 
165  float valCotTheta = cotTheta(points[0],points[1]);
166  float valEta = asinh(valCotTheta);
167  float valZip = zip(valTip, valPhi, curvature, points[0],points[1]);
168 
169  PixelTrackErrorParam param(valEta, valPt);
170  float errValPt = param.errPt();
171  float errValCot = param.errCot();
172  float errValTip = param.errTip();
173  float errValPhi = param.errPhi();
174  float errValZip = param.errZip();
175 
176 
177  float chi2 = 0;
178  if (nhits > 2) {
179  RZLine rzLine(points,errors,isBarrel);
180  float cottheta, intercept, covss, covii, covsi;
181  rzLine.fit(cottheta, intercept, covss, covii, covsi);
182  chi2 = rzLine.chi2(cottheta, intercept); //FIXME: check which intercept to use!
183  }
184 
185  PixelTrackBuilder builder;
186  Measurement1D pt(valPt, errValPt);
187  Measurement1D phi(valPhi, errValPhi);
188  Measurement1D cotTheta(valCotTheta, errValCot);
189  Measurement1D tip(valTip, errValTip);
190  Measurement1D zip(valZip, errValZip);
191 
192  return builder.build(pt, phi, cotTheta, tip, zip, chi2, iCharge, hits, theField, region.origin() );
193 }
194 
195 
196 
197 
T getParameter(std::string const &) const
edm::ESWatcher< TransientRecHitRecord > theTTRecHitBuilderWatcher
int i
Definition: DBlmapReader.cc:9
tuple cfg
Definition: looper.py:259
T perp() const
Definition: PV3DBase.h:72
GlobalPoint const & origin() const
PixelFitterByHelixProjections(const edm::ParameterSet &cfg)
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
bool isBarrel(GeomDetEnumerators::SubDetector m)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
edm::ESWatcher< TrackerDigiGeometryRecord > theTrackerWatcher
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
T y() const
Definition: PV3DBase.h:63
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
tuple zip
Definition: archive.py:476
T inversePt(T curvature, const edm::EventSetup &iSetup)
#define likely(x)
T curvature(T InversePt, const edm::EventSetup &iSetup)
edm::ESWatcher< IdealMagneticFieldRecord > theFieldWatcher
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
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
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:86
T perp2() const
Squared magnitude of transverse component.
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
const TransientTrackingRecHitBuilder * theTTRecHitBuilder
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:62
T x() const
Cartesian x coordinate.
Definition: DDAxes.h:10