CMS 3D CMS Logo

PixelFitterByHelixProjections.cc
Go to the documentation of this file.
3 
5 
9 
13 //#include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
15 
16 
19 
22 
24 
26 
32 
35 
37 
38 using namespace std;
39 
40 namespace {
41 
42  int charge(DynArray<GlobalPoint> const & points) {
43  // the cross product will tell me...
44  float dir = (points[1].x()-points[0].x())*(points[2].y()-points[1].y())
45  - (points[1].y()-points[0].y())*(points[2].x()-points[1].x());
46 
47  /*
48  GlobalVector v21 = points[1]-points[0];
49  GlobalVector v32 = points[2]-points[1];
50  float dphi = v32.phi() - v21.phi();
51  while (dphi > Geom::fpi()) dphi -= Geom::ftwoPi();
52  while (dphi < -Geom::fpi()) dphi += Geom::ftwoPi();
53  return (dphi > 0) ? -1 : 1;
54  */
55  return (dir>0) ? -1 : 1;
56  }
57 
58  float cotTheta(const GlobalPoint& inner, const GlobalPoint& outer) {
59  float dr = outer.perp()-inner.perp();
60  float dz = outer.z()-inner.z();
61  return (std::abs(dr) > 1.e-3f) ? dz/dr : 0;
62  }
63 
64  inline float phi(float xC, float yC, int charge) {
65  return (charge>0) ? std::atan2(xC,-yC) : std::atan2(-xC,yC);
66  }
67 
68  float zip(float d0, float phi_p, float curv,
69  const GlobalPoint& pinner, const GlobalPoint& pouter) {
70  //
71  //phi = asin(r*rho/2) with asin(x) ~= x+x**3/(2*3) = x(1+x*x/6);
72  //
73 
74  float phi0 = phi_p - Geom::fhalfPi();
75  GlobalPoint pca(d0*std::cos(phi0), d0*std::sin(phi0),0.);
76 
77 
78  constexpr float o24 = 1.f/24.f;
79  float rho2 = curv*curv;
80  float r1s = (pinner-pca).perp2();
81  double phi1 = std::sqrt(r1s)*(curv*0.5f)*(1.f+r1s*(rho2*o24));
82  float r2s = (pouter-pca).perp2();
83  double phi2 = std::sqrt(r2s)*(curv*0.5f)*(1.f+r2s*(rho2*o24));
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 MagneticField *field,
99  float scaleFactor):
100  theField(field),
101  thescaleErrorsForBPix1(scaleErrorsForBPix1), thescaleFactor(scaleFactor)
102 {
103  //Retrieve tracker topology from geometry
105  es->get<TrackerTopologyRcd>().get(tTopo);
106  theTopo = tTopo.product();
107 }
108 
109 std::unique_ptr<reco::Track> PixelFitterByHelixProjections::run(
110  const std::vector<const TrackingRecHit * > & hits,
111  const TrackingRegion & region,
112  const edm::EventSetup& setup) const
113 {
114  std::unique_ptr<reco::Track> ret;
115 
116  int nhits = hits.size();
117  if (nhits <2) return ret;
118 
119  declareDynArray(GlobalPoint,nhits, points);
121  declareDynArray(bool,nhits, isBarrel);
122 
123  for ( int i=0; i!=nhits; ++i) {
124  auto const & recHit = hits[i];
125  points[i] = GlobalPoint( recHit->globalPosition().basicVector()-region.origin().basicVector());
126  errors[i] = recHit->globalPositionError();
127  isBarrel[i] = recHit->detUnit()->type().isBarrel();
128  }
129 
130  CircleFromThreePoints circle = (nhits==2) ?
131  CircleFromThreePoints( GlobalPoint(0.,0.,0.), points[0], points[1]) :
132  CircleFromThreePoints(points[0],points[1],points[2]);
133 
134  float valPhi, valTip, valPt;
135 
136  int iCharge = charge(points);
137  float curvature = circle.curvature();
138 
139  if ((curvature > 1.e-4)&&
140  (LIKELY(PixelRecoUtilities::fieldInInvGev(setup)>0.01))) {
141  float invPt = PixelRecoUtilities::inversePt( circle.curvature(), setup);
142  valPt = (invPt > 1.e-4f) ? 1.f/invPt : 1.e4f;
143  CircleFromThreePoints::Vector2D center = circle.center();
144  valTip = iCharge * (center.mag()-1.f/curvature);
145  valPhi = phi(center.x(), center.y(), iCharge);
146  }
147  else {
148  valPt = 1.e4f;
149  GlobalVector direction(points[1]-points[0]);
150  valPhi = direction.barePhi();
151  valTip = -points[0].x()*sin(valPhi) + points[0].y()*cos(valPhi);
152  }
153 
154  float valCotTheta = cotTheta(points[0],points[1]);
155  float valEta = std::asinh(valCotTheta);
156  float valZip = zip(valTip, valPhi, curvature, points[0],points[1]);
157 
158  // Rescale down the error to take into accont the fact that the
159  // inner pixel barrel layer for PhaseI is closer to the interaction
160  // point. The effective scale factor has been derived by checking
161  // that the pulls of the pixelVertices derived from the pixelTracks
162  // have the correct mean and sigma.
163  float errFactor = 1.;
165  && (hits[0]->geographicalId().subdetId() == PixelSubdetector::PixelBarrel) &&
166  (theTopo->pxbLayer(hits[0]->geographicalId()) == 1))
167  errFactor = thescaleFactor;
168 
169  PixelTrackErrorParam param(valEta, valPt);
170  float errValPt = errFactor*param.errPt();
171  float errValCot = errFactor*param.errCot();
172  float errValTip = errFactor*param.errTip();
173  float errValPhi = errFactor*param.errPhi();
174  float errValZip = errFactor*param.errZip();
175 
176 
177  float chi2 = 0;
178  if (nhits > 2) {
179  RZLine rzLine(points,errors,isBarrel);
180  chi2 = rzLine.chi2();
181  }
182 
183  PixelTrackBuilder builder;
184  Measurement1D pt(valPt, errValPt);
185  Measurement1D phi(valPhi, errValPhi);
186  Measurement1D cotTheta(valCotTheta, errValCot);
187  Measurement1D tip(valTip, errValTip);
188  Measurement1D zip(valZip, errValZip);
189 
190  ret.reset(builder.build(pt, phi, cotTheta, tip, zip, chi2, iCharge, hits, theField, region.origin() ));
191  return ret;
192 }
193 
194 
195 
196 
T perp() const
Definition: PV3DBase.h:72
GlobalPoint const & origin() const
bool isBarrel(GeomDetEnumerators::SubDetector m)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:2
#define LIKELY(x)
Definition: Likely.h:20
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
T barePhi() const
Definition: PV3DBase.h:68
T inversePt(T curvature, const edm::EventSetup &iSetup)
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
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
T y() const
Cartesian y coordinate.
Definition: RZLine.h:12
unsigned int pxbLayer(const DetId &id) const
constexpr float fhalfPi()
Definition: Pi.h:37
T perp2() const
Squared magnitude of transverse component.
std::unique_ptr< reco::Track > run(const std::vector< const TrackingRecHit * > &hits, const TrackingRegion &region, const edm::EventSetup &setup) const override
T get() const
Definition: EventSetup.h:71
Definition: errors.py:1
#define declareDynArray(T, n, x)
Definition: DynArray.h:59
dbl *** dir
Definition: mlp_gen.cc:35
float fieldInInvGev(const edm::EventSetup &iSetup)
const BasicVectorType & basicVector() const
Definition: PV3DBase.h:56
float chi2() const
Definition: RZLine.h:97
T const * product() const
Definition: ESHandle.h:86
PixelFitterByHelixProjections(const edm::EventSetup *es, const MagneticField *field, bool scaleErrorsForBPix1, float scaleFactor)
#define constexpr
T x() const
Cartesian x coordinate.