CMS 3D CMS Logo

PixelFitterByHelixProjections.cc
Go to the documentation of this file.
3 
5 
9 
13 //#include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
15 
18 
21 
23 
25 
31 
33 
35 
36 using namespace std;
37 
38 namespace {
39 
40  int charge(DynArray<GlobalPoint> const& points) {
41  // the cross product will tell me...
42  float dir = (points[1].x() - points[0].x()) * (points[2].y() - points[1].y()) -
43  (points[1].y() - points[0].y()) * (points[2].x() - points[1].x());
44 
45  /*
46  GlobalVector v21 = points[1]-points[0];
47  GlobalVector v32 = points[2]-points[1];
48  float dphi = v32.phi() - v21.phi();
49  while (dphi > Geom::fpi()) dphi -= Geom::ftwoPi();
50  while (dphi < -Geom::fpi()) dphi += Geom::ftwoPi();
51  return (dphi > 0) ? -1 : 1;
52  */
53  return (dir > 0) ? -1 : 1;
54  }
55 
56  float cotTheta(const GlobalPoint& inner, const GlobalPoint& outer) {
57  float dr = outer.perp() - inner.perp();
58  float dz = outer.z() - inner.z();
59  return (std::abs(dr) > 1.e-3f) ? dz / dr : 0;
60  }
61 
62  inline float phi(float xC, float yC, int charge) { return (charge > 0) ? std::atan2(xC, -yC) : std::atan2(-xC, yC); }
63 
64  float zip(float d0, float phi_p, float curv, const GlobalPoint& pinner, const GlobalPoint& pouter) {
65  //
66  //phi = asin(r*rho/2) with asin(x) ~= x+x**3/(2*3) = x(1+x*x/6);
67  //
68 
69  float phi0 = phi_p - Geom::fhalfPi();
70  GlobalPoint pca(d0 * std::cos(phi0), d0 * std::sin(phi0), 0.);
71 
72  constexpr float o24 = 1.f / 24.f;
73  float rho2 = curv * curv;
74  float r1s = (pinner - pca).perp2();
75  double phi1 = std::sqrt(r1s) * (curv * 0.5f) * (1.f + r1s * (rho2 * o24));
76  float r2s = (pouter - pca).perp2();
77  double phi2 = std::sqrt(r2s) * (curv * 0.5f) * (1.f + r2s * (rho2 * o24));
78  double z1 = pinner.z();
79  double z2 = pouter.z();
80 
81  if (fabs(curv) > 1.e-5)
82  return z1 - phi1 / (phi1 - phi2) * (z1 - z2);
83  else {
84  double dr = std::max(std::sqrt(r2s) - std::sqrt(r1s), 1.e-5f);
85  return z1 - std::sqrt(r1s) * (z2 - z1) / dr;
86  }
87  }
88 } // namespace
89 
91  const MagneticField* field,
93  float scaleFactor)
94  : theTopo(ttopo), theField(field), thescaleErrorsForBPix1(scaleErrorsForBPix1), thescaleFactor(scaleFactor) {}
95 
96 std::unique_ptr<reco::Track> PixelFitterByHelixProjections::run(const std::vector<const TrackingRecHit*>& hits,
97  const TrackingRegion& region) const {
98  std::unique_ptr<reco::Track> ret;
99 
100  int nhits = hits.size();
101  if (nhits < 2)
102  return ret;
103 
107 
108  for (int i = 0; i != nhits; ++i) {
109  auto const& recHit = hits[i];
110  points[i] = GlobalPoint(recHit->globalPosition().basicVector() - region.origin().basicVector());
111  errors[i] = recHit->globalPositionError();
112  isBarrel[i] = recHit->detUnit()->type().isBarrel();
113  }
114 
115  CircleFromThreePoints circle = (nhits == 2) ? CircleFromThreePoints(GlobalPoint(0., 0., 0.), points[0], points[1])
116  : CircleFromThreePoints(points[0], points[1], points[2]);
117 
118  float valPhi, valTip, valPt;
119 
120  int iCharge = charge(points);
121  float curvature = circle.curvature();
122 
123  if ((curvature > 1.e-4) && (LIKELY(theField->inverseBzAtOriginInGeV()) > 0.01)) {
124  float invPt = PixelRecoUtilities::inversePt(circle.curvature(), *theField);
125  valPt = (invPt > 1.e-4f) ? 1.f / invPt : 1.e4f;
126  CircleFromThreePoints::Vector2D center = circle.center();
127  valTip = iCharge * (center.mag() - 1.f / curvature);
128  valPhi = phi(center.x(), center.y(), iCharge);
129  } else {
130  valPt = 1.e4f;
131  GlobalVector direction(points[1] - points[0]);
132  valPhi = direction.barePhi();
133  valTip = -points[0].x() * sin(valPhi) + points[0].y() * cos(valPhi);
134  }
135 
136  float valCotTheta = cotTheta(points[0], points[1]);
137  float valEta = std::asinh(valCotTheta);
138  float valZip = zip(valTip, valPhi, curvature, points[0], points[1]);
139 
140  // Rescale down the error to take into accont the fact that the
141  // inner pixel barrel layer for PhaseI is closer to the interaction
142  // point. The effective scale factor has been derived by checking
143  // that the pulls of the pixelVertices derived from the pixelTracks
144  // have the correct mean and sigma.
145  float errFactor = 1.;
146  if (thescaleErrorsForBPix1 && (hits[0]->geographicalId().subdetId() == PixelSubdetector::PixelBarrel) &&
147  (theTopo->pxbLayer(hits[0]->geographicalId()) == 1))
148  errFactor = thescaleFactor;
149 
150  PixelTrackErrorParam param(valEta, valPt);
151  float errValPt = errFactor * param.errPt();
152  float errValCot = errFactor * param.errCot();
153  float errValTip = errFactor * param.errTip();
154  float errValPhi = errFactor * param.errPhi();
155  float errValZip = errFactor * param.errZip();
156 
157  float chi2 = 0;
158  if (nhits > 2) {
159  RZLine rzLine(points, errors, isBarrel);
160  chi2 = rzLine.chi2();
161  }
162 
163  PixelTrackBuilder builder;
164  Measurement1D pt(valPt, errValPt);
165  Measurement1D phi(valPhi, errValPhi);
166  Measurement1D cotTheta(valCotTheta, errValCot);
167  Measurement1D tip(valTip, errValTip);
168  Measurement1D zip(valZip, errValZip);
169 
170  ret.reset(builder.build(pt, phi, cotTheta, tip, zip, chi2, iCharge, hits, theField, region.origin()));
171  return ret;
172 }
std::unique_ptr< reco::Track > run(const std::vector< const TrackingRecHit *> &hits, const TrackingRegion &region) const override
unsigned int pxbLayer(const DetId &id) const
T z() const
Definition: PV3DBase.h:61
ret
prodAgent to be discontinued
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float zip(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:90
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
#define LIKELY(x)
Definition: Likely.h:20
T y() const
Cartesian y coordinate.
T curvature(T InversePt, const MagneticField &field)
T perp2() const
Squared magnitude of transverse component.
T barePhi() const
Definition: PV3DBase.h:65
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
float inverseBzAtOriginInGeV() const
The inverse of field z component for this map in GeV.
Definition: MagneticField.h:52
T sqrt(T t)
Definition: SSEVec.h:23
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
constexpr G4double scaleFactor
double f[11][100]
T inversePt(T curvature, const MagneticField &field)
Definition: RZLine.h:12
constexpr float fhalfPi()
Definition: Pi.h:37
static constexpr float d0
float chi2() const
Definition: RZLine.h:95
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
float x
Definition: errors.py:1
#define declareDynArray(T, n, x)
Definition: DynArray.h:91
PixelFitterByHelixProjections(const TrackerTopology *ttopo, const MagneticField *field, bool scaleErrorsForBPix1, float scaleFactor)
T x() const
Cartesian x coordinate.