CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
KFBasedPixelFitter.cc
Go to the documentation of this file.
2 
4 
6 
9 
11 
15 
17 
20 
23 
25 
26 template <class T>
27 inline T sqr(T t) {
28  return t * t;
29 }
30 
32  : TValidTrackingRecHit(*geom) {
33  localPosition_ = LocalPoint(0., 0., 0.);
34  //neglect XY differences and BS slope
35  localError_ = LocalError(sqr(beamSpot.BeamWidthX()), 0.0, sqr(beamSpot.sigmaZ()));
36 }
37 
40  result[0] = localPosition().x();
41  return result;
42 }
44  LocalError le = localPositionError();
46  m[0][0] = le.xx();
47  return m;
48 }
50  AlgebraicMatrix matrix(1, 5, 0);
51  matrix[0][3] = 1;
52  return matrix;
53 }
54 
56  const Propagator *opropagator,
57  const TransientTrackingRecHitBuilder *ttrhBuilder,
58  const TrackerGeometry *tracker,
59  const MagneticField *field,
60  const reco::BeamSpot *beamSpot)
61  : thePropagator(propagator),
62  theOPropagator(opropagator),
63  theTTRHBuilder(ttrhBuilder),
64  theTracker(tracker),
65  theField(field),
66  theBeamSpot(beamSpot) {}
67 
68 std::unique_ptr<reco::Track> KFBasedPixelFitter::run(const std::vector<const TrackingRecHit *> &hits,
69  const TrackingRegion &region) const {
70  std::unique_ptr<reco::Track> ret;
71 
72  int nhits = hits.size();
73  if (nhits < 2)
74  return ret;
75 
76  float ptMin = region.ptMin();
77 
78  const GlobalPoint &vertexPos = region.origin();
79  GlobalError vertexErr(sqr(region.originRBound()), 0, sqr(region.originRBound()), 0, 0, sqr(region.originZBound()));
80 
81  std::vector<GlobalPoint> points(nhits);
82  points[0] = theTracker->idToDet(hits[0]->geographicalId())->toGlobal(hits[0]->localPosition());
83  points[1] = theTracker->idToDet(hits[1]->geographicalId())->toGlobal(hits[1]->localPosition());
84  points[2] = theTracker->idToDet(hits[2]->geographicalId())->toGlobal(hits[2]->localPosition());
85 
86  //
87  //initial Kinematics
88  //
89  GlobalVector initMom;
90  int charge;
91  float theta;
92  CircleFromThreePoints circle(points[0], points[1], points[2]);
93  if (circle.curvature() > 1.e-4) {
94  float invPt = PixelRecoUtilities::inversePt(circle.curvature(), *theField);
95  float valPt = 1.f / invPt;
96  float chargeTmp = (points[1].x() - points[0].x()) * (points[2].y() - points[1].y()) -
97  (points[1].y() - points[0].y()) * (points[2].x() - points[1].x());
98  charge = (chargeTmp > 0) ? -1 : 1;
99  float valPhi = (charge > 0) ? std::atan2(circle.center().x(), -circle.center().y())
100  : std::atan2(-circle.center().x(), circle.center().y());
101  theta = GlobalVector(points[1] - points[0]).theta();
102  initMom = GlobalVector(valPt * cos(valPhi), valPt * sin(valPhi), valPt / tan(theta));
103  } else {
104  initMom = GlobalVector(points[1] - points[0]);
105  initMom *= 10000. / initMom.perp();
106  charge = 1;
107  theta = initMom.theta();
108  }
109  GlobalTrajectoryParameters initialKine(vertexPos, initMom, TrackCharge(charge), theField);
110 
111  //
112  // initial error
113  //
114  AlgebraicSymMatrix55 C = ROOT::Math::SMatrixIdentity();
115  float sin2th = sqr(sin(theta));
116  float minC00 = 1.0;
117  C[0][0] = std::max(sin2th / sqr(ptMin), minC00);
118  float zErr = vertexErr.czz();
119  float transverseErr = vertexErr.cxx(); // assume equal cxx cyy
120  C[3][3] = transverseErr;
121  C[4][4] = zErr * sin2th + transverseErr * (1 - sin2th);
122  CurvilinearTrajectoryError initialError(C);
123 
124  FreeTrajectoryState fts(initialKine, initialError);
125 
126  // get updator
128 
129  // Now update initial state track using information from hits.
130  TrajectoryStateOnSurface outerState;
131  DetId outerDetId = 0;
132  const TrackingRecHit *hit = nullptr;
133  for (unsigned int iHit = 0; iHit < hits.size(); iHit++) {
134  hit = hits[iHit];
135  if (iHit == 0)
136  outerState = thePropagator->propagate(fts, theTracker->idToDet(hit->geographicalId())->surface());
137  outerDetId = hit->geographicalId();
139  if (!state.isValid())
140  return ret;
141  // TransientTrackingRecHit::RecHitPointer recHit = (theTTRHBuilder->build(hit))->clone(state);
143  outerState = updator.update(state, *recHit);
144  if (!outerState.isValid())
145  return ret;
146  }
147 
148  TrajectoryStateOnSurface innerState = outerState;
149  DetId innerDetId = 0;
150  innerState.rescaleError(100000.);
151  for (int iHit = 2; iHit >= 0; --iHit) {
152  hit = hits[iHit];
153  innerDetId = hit->geographicalId();
155  if (!state.isValid())
156  return ret;
157  // TransientTrackingRecHit::RecHitPointer recHit = (theTTRHBuilder->build(hit))->clone(state);
159  innerState = updator.update(state, *recHit);
160  if (!innerState.isValid())
161  return ret;
162  }
163 
164  // extrapolate to vertex
165  auto impactPointState = TransverseImpactPointExtrapolator(theField).extrapolate(innerState, vertexPos);
166  if (!impactPointState.isValid())
167  return ret;
168 
169  //
170  // optionally update impact point state with Bs constraint
171  // using this potion makes sense if vertexPos (from TrackingRegion is centerewd at BeamSpot).
172  //
173  if (theBeamSpot) {
174  MyBeamSpotGeomDet bsgd(Plane::build(impactPointState.surface().position(), impactPointState.surface().rotation()));
175  MyBeamSpotHit bsrh(*theBeamSpot, &bsgd);
176  impactPointState = updator.update(impactPointState, bsrh); //update
177  impactPointState =
178  TransverseImpactPointExtrapolator(theField).extrapolate(impactPointState, vertexPos); //reextrapolate
179  if (!impactPointState.isValid())
180  return ret;
181  }
182 
183  int ndof = 2 * hits.size() - 5;
184  GlobalPoint vv = impactPointState.globalPosition();
185  math::XYZPoint pos(vv.x(), vv.y(), vv.z());
186  GlobalVector pp = impactPointState.globalMomentum();
187  math::XYZVector mom(pp.x(), pp.y(), pp.z());
188 
189  float chi2 = 0.;
190  ret = std::make_unique<reco::Track>(
191  chi2, ndof, pos, mom, impactPointState.charge(), impactPointState.curvilinearError());
192 
193  /*
194  vv = outerState.globalPosition();
195  pp = outerState.globalMomentum();
196  math::XYZPoint outerPosition( vv.x(), vv.y(), vv.z());
197  math::XYZVector outerMomentum( pp.x(), pp.y(), pp.z());
198  vv = innerState.globalPosition();
199  pp = innerState.globalMomentum();
200  math::XYZPoint innerPosition( vv.x(), vv.y(), vv.z());
201  math::XYZVector innerMomentum( pp.x(), pp.y(), pp.z());
202 
203  reco::TrackExtra extra( outerPosition, outerMomentum, true,
204  innerPosition, innerMomentum, true,
205  outerState.curvilinearError(), outerDetId,
206  innerState.curvilinearError(), innerDetId,
207  anyDirection);
208 */
209 
210  // std::cout <<"TRACK CREATED" << std::endl;
211  return ret;
212 }
float originRBound() const
bounds the particle vertex in the transverse plane
const Propagator * theOPropagator
tuple propagator
float xx() const
Definition: LocalError.h:22
tuple ret
prodAgent to be discontinued
const Propagator * thePropagator
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
T perp() const
Definition: PV3DBase.h:69
GlobalPoint const & origin() const
tuple pp
Definition: createTree.py:17
AlgebraicSymMatrix parametersError() const override
const TrackerGeometry * theTracker
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
constexpr float ptMin
Geom::Theta< T > theta() const
T y() const
Definition: PV3DBase.h:60
int sqr(const T &t)
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
const MagneticField * theField
tuple result
Definition: mps_fire.py:311
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72
int TrackCharge
Definition: TrackCharge.h:4
AlgebraicMatrix projectionMatrix() const override
TrajectoryStateOnSurface update(const TrajectoryStateOnSurface &, const TrackingRecHit &) const override
Definition: KFUpdator.cc:177
CLHEP::HepMatrix AlgebraicMatrix
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
static PlanePointer build(Args &&...args)
Definition: Plane.h:33
T z() const
Definition: PV3DBase.h:61
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
MyBeamSpotHit(const reco::BeamSpot &beamSpot, const GeomDet *geom)
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
double BeamWidthX() const
beam width X
Definition: BeamSpot.h:82
std::unique_ptr< reco::Track > run(const std::vector< const TrackingRecHit * > &hits, const TrackingRegion &region) const override
T inversePt(T curvature, const MagneticField &field)
T y() const
Cartesian y coordinate.
AlgebraicVector parameters() const override
const TrackerGeomDet * idToDet(DetId) const override
float originZBound() const
bounds the particle vertex in the longitudinal plane
std::shared_ptr< TrackingRecHit const > RecHitPointer
const reco::BeamSpot * theBeamSpot
KFBasedPixelFitter(const Propagator *propagator, const Propagator *opropagator, const TransientTrackingRecHitBuilder *ttrhBuilder, const TrackerGeometry *tracker, const MagneticField *field, const reco::BeamSpot *beamSpot)
const TransientTrackingRecHitBuilder * theTTRHBuilder
Definition: DetId.h:17
CLHEP::HepVector AlgebraicVector
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
double sigmaZ() const
sigma z
Definition: BeamSpot.h:76
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
float ptMin() const
minimal pt of interest
TrajectoryStateOnSurface propagate(STA const &state, SUR const &surface) const
Definition: Propagator.h:50
TrajectoryStateOnSurface extrapolate(const FreeTrajectoryState &fts, const GlobalPoint &vtx) const
extrapolation with default (=geometrical) propagator
CLHEP::HepSymMatrix AlgebraicSymMatrix
DetId geographicalId() const
long double T
T x() const
Definition: PV3DBase.h:59
T x() const
Cartesian x coordinate.
Global3DVector GlobalVector
Definition: GlobalVector.h:10