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"
39 
40 using namespace std;
41 
42 namespace {
43 
44  int charge(const std::vector<GlobalPoint> & points) {
45  // the cross product will tell me...
46  float dir = (points[1].x()-points[0].x())*(points[2].y()-points[1].y())
47  - (points[1].y()-points[0].y())*(points[2].x()-points[1].x());
48 
49  /*
50  GlobalVector v21 = points[1]-points[0];
51  GlobalVector v32 = points[2]-points[1];
52  float dphi = v32.phi() - v21.phi();
53  while (dphi > Geom::fpi()) dphi -= Geom::ftwoPi();
54  while (dphi < -Geom::fpi()) dphi += Geom::ftwoPi();
55  return (dphi > 0) ? -1 : 1;
56  */
57  return (dir>0) ? -1 : 1;
58  }
59 
60  float cotTheta(const GlobalPoint& inner, const GlobalPoint& outer) {
61  float dr = outer.perp()-inner.perp();
62  float dz = outer.z()-inner.z();
63  return (std::abs(dr) > 1.e-3f) ? dz/dr : 0;
64  }
65 
66  inline float phi(float xC, float yC, int charge) {
67  return (charge>0) ? std::atan2(xC,-yC) : std::atan2(-xC,yC);
68  }
69 
70  float zip(float d0, float phi_p, float curv,
71  const GlobalPoint& pinner, const GlobalPoint& pouter) {
72  //
73  //phi = asin(r*rho/2) with asin(x) ~= x+x**3/(2*3) = x(1+x*x/6);
74  //
75 
76  float phi0 = phi_p - Geom::fhalfPi();
77  GlobalPoint pca(d0*std::cos(phi0), d0*std::sin(phi0),0.);
78 
79  float rho2 = curv*curv;
80  float r1s = (pinner-pca).perp2();
81  double phi1 = std::sqrt(r1s)*(curv*0.5f)*(1.f+r1s*(rho2/24.f));
82  float r2s = (pouter-pca).perp2();
83  double phi2 = std::sqrt(r2s)*(curv*0.5f)*(1.f+r2s*(rho2/24.f));
84  double z1 = pinner.z();
85  double z2 = pouter.z();
86 
87  if (fabs(curv)>1.e-5)
88  return z1 - phi1/(phi1-phi2)*(z1-z2);
89  else {
90  double dr = std::max(std::sqrt(r2s)-std::sqrt(r1s),1.e-5f);
91  return z1-std::sqrt(r1s)*(z2-z1)/dr;
92  }
93  }
94 }
95 
97  const edm::ParameterSet& cfg)
98  : theConfig(cfg), theTracker(0), theField(0), theTTRecHitBuilder(0) {}
99 
101  const edm::EventSetup& es,
102  const std::vector<const TrackingRecHit * > & hits,
103  const TrackingRegion & region) const
104 {
105  int nhits = hits.size();
106  if (nhits <2) return 0;
107 
108  vector<GlobalPoint> points(nhits);
109  vector<GlobalError> errors(nhits);
110  vector<bool> isBarrel(nhits);
111 
112  static edm::ESWatcher<TrackerDigiGeometryRecord> watcherTrackerDigiGeometryRecord;
113  if (!theTracker || watcherTrackerDigiGeometryRecord.check(es)) {
115  es.get<TrackerDigiGeometryRecord>().get(trackerESH);
116  theTracker = trackerESH.product();
117  }
118 
119  static edm::ESWatcher<IdealMagneticFieldRecord> watcherIdealMagneticFieldRecord;
120  if (!theField || watcherIdealMagneticFieldRecord.check(es)) {
122  es.get<IdealMagneticFieldRecord>().get(fieldESH);
123  theField = fieldESH.product();
124  }
125 
126  static edm::ESWatcher<TransientRecHitRecord> watcherTransientRecHitRecord;
127  if (!theTTRecHitBuilder || watcherTransientRecHitRecord.check(es)) {
129  std::string builderName = theConfig.getParameter<std::string>("TTRHBuilder");
130  es.get<TransientRecHitRecord>().get(builderName,ttrhbESH);
131  theTTRecHitBuilder = ttrhbESH.product();
132  }
133 
134 
135  for ( int i=0; i!=nhits; ++i) {
137  points[i] = GlobalPoint( recHit->globalPosition().x()-region.origin().x(),
138  recHit->globalPosition().y()-region.origin().y(),
139  recHit->globalPosition().z()-region.origin().z()
140  );
141  errors[i] = recHit->globalPositionError();
142  isBarrel[i] = recHit->detUnit()->type().isBarrel();
143  }
144 
145  CircleFromThreePoints circle = (nhits==2) ?
146  CircleFromThreePoints( GlobalPoint(0.,0.,0.), points[0], points[1]) :
147  CircleFromThreePoints(points[0],points[1],points[2]);
148 
149  float valPhi, valTip, valPt;
150 
151  int iCharge = charge(points);
152  float curvature = circle.curvature();
153 
154  if ((curvature > 1.e-4)&&
155  (likely(theField->inTesla(GlobalPoint(0.,0.,0.)).z()>0.01))) {
156  float invPt = PixelRecoUtilities::inversePt( circle.curvature(), es);
157  valPt = (invPt > 1.e-4f) ? 1.f/invPt : 1.e4f;
158  CircleFromThreePoints::Vector2D center = circle.center();
159  valTip = iCharge * (center.mag()-1.f/curvature);
160  valPhi = phi(center.x(), center.y(), iCharge);
161  }
162  else {
163  valPt = 1.e4f;
164  GlobalVector direction(points[1]-points[0]);
165  valPhi = direction.phi();
166  valTip = -points[0].x()*sin(valPhi) + points[0].y()*cos(valPhi);
167  }
168 
169  float valCotTheta = cotTheta(points[0],points[1]);
170  float valEta = asinh(valCotTheta);
171  float valZip = zip(valTip, valPhi, curvature, points[0],points[1]);
172 
173  PixelTrackErrorParam param(valEta, valPt);
174  float errValPt = param.errPt();
175  float errValCot = param.errCot();
176  float errValTip = param.errTip();
177  float errValPhi = param.errPhi();
178  float errValZip = param.errZip();
179 
180 
181  float chi2 = 0;
182  if (nhits > 2) {
183  RZLine rzLine(points,errors,isBarrel);
184  float cottheta, intercept, covss, covii, covsi;
185  rzLine.fit(cottheta, intercept, covss, covii, covsi);
186  chi2 = rzLine.chi2(cottheta, intercept); //FIXME: check which intercept to use!
187  }
188 
189  PixelTrackBuilder builder;
190  Measurement1D pt(valPt, errValPt);
191  Measurement1D phi(valPhi, errValPhi);
192  Measurement1D cotTheta(valCotTheta, errValCot);
193  Measurement1D tip(valTip, errValTip);
194  Measurement1D zip(valZip, errValZip);
195 
196  return builder.build(pt, phi, cotTheta, tip, zip, chi2, iCharge, hits, theField, region.origin() );
197 }
198 
199 
200 
201 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
T perp() const
Definition: PV3DBase.h:71
virtual GlobalPoint origin() const =0
PixelFitterByHelixProjections(const edm::ParameterSet &cfg)
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Phi< T > phi() const
Definition: PV3DBase.h:68
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
T y() const
Definition: PV3DBase.h:62
#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
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:46
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
T z() const
Definition: PV3DBase.h:63
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]
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
#define likely(x)
Definition: Likely.h:20
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:61
T x() const
Cartesian x coordinate.
Definition: DDAxes.h:10