CMS 3D CMS Logo

KFBasedPixelFitter.cc
Go to the documentation of this file.
2 
6 
8 
11 
13 
17 
20 
23 
27 
31 
32 /*
33 #include "RecoVertex/KalmanVertexFit/interface/SingleTrackVertexConstraint.h"
34 #include "TrackingTools/TransientTrack/interface/TransientTrackFromFTSFactory.h"
35 
36 #include "TrackingTools/TransientTrack/interface/TransientTrackFromFTSFactory.h"
37 #include "Alignment/ReferenceTrajectories/interface/BeamSpotTransientTrackingRecHit.h"
38 #include "Alignment/ReferenceTrajectories/interface/BeamSpotGeomDet.h"
39 #include <TrackingTools/PatternTools/interface/TSCPBuilderNoMaterial.h>
40 #include <TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h>
41 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateClosestToPoint.h"
42 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateClosestToBeamLine.h"
43 */
44 
45 #include <sstream>
46 template <class T>
47 inline T sqr(T t) {
48  return t * t;
49 }
50 
52  : TValidTrackingRecHit(*geom) {
53  localPosition_ = LocalPoint(0., 0., 0.);
54  //neglect XY differences and BS slope
55  localError_ = LocalError(sqr(beamSpot.BeamWidthX()), 0.0, sqr(beamSpot.sigmaZ()));
56 }
57 
60  result[0] = localPosition().x();
61  return result;
62 }
66  m[0][0] = le.xx();
67  return m;
68 }
70  AlgebraicMatrix matrix(1, 5, 0);
71  matrix[0][3] = 1;
72  return matrix;
73 }
74 
76  const Propagator *opropagator,
77  const TransientTrackingRecHitBuilder *ttrhBuilder,
78  const TrackerGeometry *tracker,
79  const MagneticField *field,
80  const reco::BeamSpot *beamSpot)
81  : thePropagator(propagator),
82  theOPropagator(opropagator),
83  theTTRHBuilder(ttrhBuilder),
84  theTracker(tracker),
85  theField(field),
86  theBeamSpot(beamSpot) {}
87 
88 std::unique_ptr<reco::Track> KFBasedPixelFitter::run(const std::vector<const TrackingRecHit *> &hits,
89  const TrackingRegion &region,
90  const edm::EventSetup &setup) const {
91  std::unique_ptr<reco::Track> ret;
92 
93  int nhits = hits.size();
94  if (nhits < 2)
95  return ret;
96 
97  float ptMin = region.ptMin();
98 
99  const GlobalPoint &vertexPos = region.origin();
100  GlobalError vertexErr(sqr(region.originRBound()), 0, sqr(region.originRBound()), 0, 0, sqr(region.originZBound()));
101 
102  std::vector<GlobalPoint> points(nhits);
103  points[0] = theTracker->idToDet(hits[0]->geographicalId())->toGlobal(hits[0]->localPosition());
104  points[1] = theTracker->idToDet(hits[1]->geographicalId())->toGlobal(hits[1]->localPosition());
105  points[2] = theTracker->idToDet(hits[2]->geographicalId())->toGlobal(hits[2]->localPosition());
106 
107  //
108  //initial Kinematics
109  //
110  GlobalVector initMom;
111  int charge;
112  float theta;
113  CircleFromThreePoints circle(points[0], points[1], points[2]);
114  if (circle.curvature() > 1.e-4) {
115  float invPt = PixelRecoUtilities::inversePt(circle.curvature(), setup);
116  float valPt = 1.f / invPt;
117  float chargeTmp = (points[1].x() - points[0].x()) * (points[2].y() - points[1].y()) -
118  (points[1].y() - points[0].y()) * (points[2].x() - points[1].x());
119  charge = (chargeTmp > 0) ? -1 : 1;
120  float valPhi = (charge > 0) ? std::atan2(circle.center().x(), -circle.center().y())
121  : std::atan2(-circle.center().x(), circle.center().y());
122  theta = GlobalVector(points[1] - points[0]).theta();
123  initMom = GlobalVector(valPt * cos(valPhi), valPt * sin(valPhi), valPt / tan(theta));
124  } else {
125  initMom = GlobalVector(points[1] - points[0]);
126  initMom *= 10000. / initMom.perp();
127  charge = 1;
128  theta = initMom.theta();
129  }
130  GlobalTrajectoryParameters initialKine(vertexPos, initMom, TrackCharge(charge), theField);
131 
132  //
133  // initial error
134  //
135  AlgebraicSymMatrix55 C = ROOT::Math::SMatrixIdentity();
136  float sin2th = sqr(sin(theta));
137  float minC00 = 1.0;
138  C[0][0] = std::max(sin2th / sqr(ptMin), minC00);
139  float zErr = vertexErr.czz();
140  float transverseErr = vertexErr.cxx(); // assume equal cxx cyy
141  C[3][3] = transverseErr;
142  C[4][4] = zErr * sin2th + transverseErr * (1 - sin2th);
143  CurvilinearTrajectoryError initialError(C);
144 
145  FreeTrajectoryState fts(initialKine, initialError);
146 
147  // get updator
149 
150  // Now update initial state track using information from hits.
151  TrajectoryStateOnSurface outerState;
152  DetId outerDetId = 0;
153  const TrackingRecHit *hit = nullptr;
154  for (unsigned int iHit = 0; iHit < hits.size(); iHit++) {
155  hit = hits[iHit];
156  if (iHit == 0)
157  outerState = thePropagator->propagate(fts, theTracker->idToDet(hit->geographicalId())->surface());
158  outerDetId = hit->geographicalId();
159  TrajectoryStateOnSurface state = thePropagator->propagate(outerState, theTracker->idToDet(outerDetId)->surface());
160  if (!state.isValid())
161  return ret;
162  // TransientTrackingRecHit::RecHitPointer recHit = (theTTRHBuilder->build(hit))->clone(state);
164  outerState = updator.update(state, *recHit);
165  if (!outerState.isValid())
166  return ret;
167  }
168 
169  TrajectoryStateOnSurface innerState = outerState;
170  DetId innerDetId = 0;
171  innerState.rescaleError(100000.);
172  for (int iHit = 2; iHit >= 0; --iHit) {
173  hit = hits[iHit];
174  innerDetId = hit->geographicalId();
175  TrajectoryStateOnSurface state = theOPropagator->propagate(innerState, theTracker->idToDet(innerDetId)->surface());
176  if (!state.isValid())
177  return ret;
178  // TransientTrackingRecHit::RecHitPointer recHit = (theTTRHBuilder->build(hit))->clone(state);
180  innerState = updator.update(state, *recHit);
181  if (!innerState.isValid())
182  return ret;
183  }
184 
185  // extrapolate to vertex
186  TrajectoryStateOnSurface impactPointState =
188  if (!impactPointState.isValid())
189  return ret;
190 
191  //
192  // optionally update impact point state with Bs constraint
193  // using this potion makes sense if vertexPos (from TrackingRegion is centerewd at BeamSpot).
194  //
195  if (theBeamSpot) {
196  MyBeamSpotGeomDet bsgd(Plane::build(impactPointState.surface().position(), impactPointState.surface().rotation()));
197  MyBeamSpotHit bsrh(*theBeamSpot, &bsgd);
198  impactPointState = updator.update(impactPointState, bsrh); //update
199  impactPointState =
200  TransverseImpactPointExtrapolator(theField).extrapolate(impactPointState, vertexPos); //reextrapolate
201  if (!impactPointState.isValid())
202  return ret;
203  }
204 
205  int ndof = 2 * hits.size() - 5;
206  GlobalPoint vv = impactPointState.globalPosition();
207  math::XYZPoint pos(vv.x(), vv.y(), vv.z());
208  GlobalVector pp = impactPointState.globalMomentum();
209  math::XYZVector mom(pp.x(), pp.y(), pp.z());
210 
211  float chi2 = 0.;
212  ret = std::make_unique<reco::Track>(
213  chi2, ndof, pos, mom, impactPointState.charge(), impactPointState.curvilinearError());
214 
215  /*
216  vv = outerState.globalPosition();
217  pp = outerState.globalMomentum();
218  math::XYZPoint outerPosition( vv.x(), vv.y(), vv.z());
219  math::XYZVector outerMomentum( pp.x(), pp.y(), pp.z());
220  vv = innerState.globalPosition();
221  pp = innerState.globalMomentum();
222  math::XYZPoint innerPosition( vv.x(), vv.y(), vv.z());
223  math::XYZVector innerMomentum( pp.x(), pp.y(), pp.z());
224 
225  reco::TrackExtra extra( outerPosition, outerMomentum, true,
226  innerPosition, innerMomentum, true,
227  outerState.curvilinearError(), outerDetId,
228  innerState.curvilinearError(), innerDetId,
229  anyDirection);
230 */
231 
232  // std::cout <<"TRACK CREATED" << std::endl;
233  return ret;
234 }
float originRBound() const
bounds the particle vertex in the transverse plane
const Propagator * theOPropagator
float xx() const
Definition: LocalError.h:22
const Propagator * thePropagator
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
T perp() const
Definition: PV3DBase.h:69
GlobalPoint const & origin() const
std::unique_ptr< reco::Track > run(const std::vector< const TrackingRecHit * > &hits, const TrackingRegion &region, const edm::EventSetup &setup) const override
ret
prodAgent to be discontinued
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
const CurvilinearTrajectoryError & curvilinearError() const
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
T y() const
Definition: PV3DBase.h:60
GlobalPoint globalPosition() const
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
const MagneticField * theField
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72
T inversePt(T curvature, const edm::EventSetup &iSetup)
int TrackCharge
Definition: TrackCharge.h:4
TrajectoryStateOnSurface update(const TrajectoryStateOnSurface &, const TrackingRecHit &) const override
Definition: KFUpdator.cc:177
const SurfaceType & surface() const
T sqr(T t)
CLHEP::HepMatrix AlgebraicMatrix
AlgebraicVector parameters() const override
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
T y() const
Cartesian y coordinate.
float originZBound() const
bounds the particle vertex in the longitudinal plane
std::shared_ptr< TrackingRecHit const > RecHitPointer
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
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
AlgebraicMatrix projectionMatrix() const override
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
GlobalVector globalMomentum() const
CLHEP::HepSymMatrix AlgebraicSymMatrix
const TrackerGeomDet * idToDet(DetId) const override
AlgebraicSymMatrix parametersError() const override
const RotationType & rotation() const
DetId geographicalId() const
LocalError localPositionError() const override
long double T
T x() const
Definition: PV3DBase.h:59
const PositionType & position() const
T x() const
Cartesian x coordinate.
Global3DVector GlobalVector
Definition: GlobalVector.h:10
LocalPoint localPosition() const override