CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Private Attributes
KFBasedPixelFitter Class Reference

#include <KFBasedPixelFitter.h>

Inheritance diagram for KFBasedPixelFitter:
PixelFitterBase

Classes

class  MyBeamSpotGeomDet
 
class  MyBeamSpotHit
 

Public Member Functions

 KFBasedPixelFitter (const edm::EventSetup *es, const Propagator *propagator, const Propagator *opropagator, const TransientTrackingRecHitBuilder *ttrhBuilder, const TrackerGeometry *tracker, const MagneticField *field, const reco::BeamSpot *beamSpot)
 
std::unique_ptr< reco::Trackrun (const std::vector< const TrackingRecHit * > &hits, const TrackingRegion &region) const override
 
 ~KFBasedPixelFitter () override
 
- Public Member Functions inherited from PixelFitterBase
virtual reco::Trackrun (const edm::EventSetup &es, const std::vector< const TrackingRecHit * > &hits, const TrackingRegion &region) const
 
virtual reco::Trackrun (const edm::Event &ev, const edm::EventSetup &es, const std::vector< const TrackingRecHit * > &hits, const TrackingRegion &region) const
 
virtual ~PixelFitterBase ()
 

Private Attributes

const reco::BeamSpottheBeamSpot
 
const edm::EventSetuptheES
 
const MagneticFieldtheField
 
const PropagatortheOPropagator
 
const PropagatorthePropagator
 
const TrackerGeometrytheTracker
 
const TransientTrackingRecHitBuildertheTTRHBuilder
 

Detailed Description

Definition at line 22 of file KFBasedPixelFitter.h.

Constructor & Destructor Documentation

KFBasedPixelFitter::KFBasedPixelFitter ( const edm::EventSetup es,
const Propagator propagator,
const Propagator opropagator,
const TransientTrackingRecHitBuilder ttrhBuilder,
const TrackerGeometry tracker,
const MagneticField field,
const reco::BeamSpot beamSpot 
)

Definition at line 78 of file KFBasedPixelFitter.cc.

81  :
82  theES(es), thePropagator(propagator), theOPropagator(opropagator), theTTRHBuilder(ttrhBuilder),
83  theTracker(tracker), theField(field), theBeamSpot(beamSpot)
84 {}
const Propagator * theOPropagator
const Propagator * thePropagator
const TrackerGeometry * theTracker
const MagneticField * theField
const reco::BeamSpot * theBeamSpot
const edm::EventSetup * theES
const TransientTrackingRecHitBuilder * theTTRHBuilder
KFBasedPixelFitter::~KFBasedPixelFitter ( )
inlineoverride

Definition at line 28 of file KFBasedPixelFitter.h.

References hfClusterShapes_cfi::hits, and findQualityFiles::run.

28 {}

Member Function Documentation

std::unique_ptr< reco::Track > KFBasedPixelFitter::run ( const std::vector< const TrackingRecHit * > &  hits,
const TrackingRegion region 
) const
overridevirtual

Reimplemented from PixelFitterBase.

Definition at line 86 of file KFBasedPixelFitter.cc.

References TransientTrackingRecHitBuilder::build(), Plane::build(), patCaloMETCorrections_cff::C, CircleFromThreePoints::center(), ALCARECOTkAlJpsiMuMu_cff::charge, TrajectoryStateOnSurface::charge(), vertices_cff::chi2, funct::cos(), CircleFromThreePoints::curvature(), TrajectoryStateOnSurface::curvilinearError(), GlobalErrorBase< T, ErrorWeightType >::cxx(), GlobalErrorBase< T, ErrorWeightType >::czz(), TransverseImpactPointExtrapolator::extrapolate(), TrackingRecHit::geographicalId(), TrajectoryStateOnSurface::globalMomentum(), TrajectoryStateOnSurface::globalPosition(), TrackerGeometry::idToDet(), PixelRecoUtilities::inversePt(), TrajectoryStateOnSurface::isValid(), SiStripPI::max, ndof, nhits, TrackingRegion::origin(), TrackingRegion::originRBound(), TrackingRegion::originZBound(), PV3DBase< T, PVType, FrameType >::perp(), hiPixelPairStep_cff::points, GloballyPositioned< T >::position(), createTree::pp, Propagator::propagate(), TrackingRegion::ptMin(), ptMin, rpcPointValidation_cfi::recHit, TrajectoryStateOnSurface::rescaleError(), GloballyPositioned< T >::rotation(), funct::sin(), sqr(), GeomDet::surface(), TrajectoryStateOnSurface::surface(), funct::tan(), theBeamSpot, theES, theField, theOPropagator, thePropagator, PV3DBase< T, PVType, FrameType >::theta(), theta(), theTracker, theTTRHBuilder, GeomDet::toGlobal(), KFUpdator::update(), gsfElectronCkfTrackCandidateMaker_cff::updator, x, PV3DBase< T, PVType, FrameType >::x(), Basic2DVector< T >::x(), y, PV3DBase< T, PVType, FrameType >::y(), Basic2DVector< T >::y(), and PV3DBase< T, PVType, FrameType >::z().

86  {
87  std::unique_ptr<reco::Track> ret;
88 
89  int nhits = hits.size();
90  if (nhits <2) return ret;
91 
92 
93  float ptMin = region.ptMin();
94 
95  const GlobalPoint& vertexPos = region.origin();
96  GlobalError vertexErr( sqr(region.originRBound()), 0, sqr(region.originRBound()), 0, 0, sqr(region.originZBound()));
97 
98  std::vector<GlobalPoint> points(nhits);
99  points[0] = theTracker->idToDet(hits[0]->geographicalId())->toGlobal(hits[0]->localPosition());
100  points[1] = theTracker->idToDet(hits[1]->geographicalId())->toGlobal(hits[1]->localPosition());
101  points[2] = theTracker->idToDet(hits[2]->geographicalId())->toGlobal(hits[2]->localPosition());
102 
103  //
104  //initial Kinematics
105  //
106  GlobalVector initMom;
107  int charge;
108  float theta;
109  CircleFromThreePoints circle(points[0], points[1], points[2]);
110  if (circle.curvature() > 1.e-4) {
111  float invPt = PixelRecoUtilities::inversePt( circle.curvature(), *theES);
112  float valPt = 1.f/invPt;
113  float chargeTmp = (points[1].x()-points[0].x())*(points[2].y()-points[1].y())
114  - (points[1].y()-points[0].y())*(points[2].x()-points[1].x());
115  charge = (chargeTmp>0) ? -1 : 1;
116  float valPhi = (charge>0) ? std::atan2(circle.center().x(),-circle.center().y()) : std::atan2(-circle.center().x(),circle.center().y());
117  theta = GlobalVector(points[1]-points[0]).theta();
118  initMom = GlobalVector(valPt*cos(valPhi), valPt*sin(valPhi), valPt/tan(theta));
119  }
120  else {
121  initMom = GlobalVector(points[1]-points[0]);
122  initMom *= 10000./initMom.perp();
123  charge = 1;
124  theta = initMom.theta();
125  }
126  GlobalTrajectoryParameters initialKine(vertexPos, initMom, TrackCharge(charge), theField);
127 
128  //
129  // initial error
130  //
131  AlgebraicSymMatrix55 C = ROOT::Math::SMatrixIdentity();
132  float sin2th = sqr(sin(theta));
133  float minC00 = 1.0;
134  C[0][0] = std::max(sin2th/sqr(ptMin), minC00);
135  float zErr = vertexErr.czz();
136  float transverseErr = vertexErr.cxx(); // assume equal cxx cyy
137  C[3][3] = transverseErr;
138  C[4][4] = zErr*sin2th + transverseErr*(1-sin2th);
139  CurvilinearTrajectoryError initialError(C);
140 
141  FreeTrajectoryState fts(initialKine, initialError);
142 
143  // get updator
145 
146  // Now update initial state track using information from hits.
147  TrajectoryStateOnSurface outerState;
148  DetId outerDetId = 0;
149  const TrackingRecHit* hit = nullptr;
150  for ( unsigned int iHit = 0; iHit < hits.size(); iHit++) {
151  hit = hits[iHit];
152  if (iHit==0) outerState = thePropagator->propagate(fts,theTracker->idToDet(hit->geographicalId())->surface());
153  outerDetId = hit->geographicalId();
154  TrajectoryStateOnSurface state = thePropagator->propagate(outerState, theTracker->idToDet(outerDetId)->surface());
155  if (!state.isValid()) return ret;
156 // TransientTrackingRecHit::RecHitPointer recHit = (theTTRHBuilder->build(hit))->clone(state);
158  outerState = updator.update(state, *recHit);
159  if (!outerState.isValid()) return ret;
160  }
161 
162 
163 
164  TrajectoryStateOnSurface innerState = outerState;
165  DetId innerDetId = 0;
166  innerState.rescaleError(100000.);
167  for ( int iHit = 2; iHit >= 0; --iHit) {
168  hit = hits[iHit];
169  innerDetId = hit->geographicalId();
170  TrajectoryStateOnSurface state = theOPropagator->propagate(innerState, theTracker->idToDet(innerDetId)->surface());
171  if (!state.isValid()) return ret;
172 // TransientTrackingRecHit::RecHitPointer recHit = (theTTRHBuilder->build(hit))->clone(state);
174  innerState = updator.update(state, *recHit);
175  if (!innerState.isValid()) return ret;
176  }
177 
178 
179  // extrapolate to vertex
180  TrajectoryStateOnSurface impactPointState = TransverseImpactPointExtrapolator(theField).extrapolate( innerState, vertexPos);
181  if (!impactPointState.isValid()) return ret;
182 
183  //
184  // optionally update impact point state with Bs constraint
185  // using this potion makes sense if vertexPos (from TrackingRegion is centerewd at BeamSpot).
186  //
187  if (theBeamSpot) {
188  MyBeamSpotGeomDet bsgd(Plane::build(impactPointState.surface().position(), impactPointState.surface().rotation()));
189  MyBeamSpotHit bsrh(*theBeamSpot, &bsgd);
190  impactPointState = updator.update(impactPointState, bsrh); //update
191  impactPointState = TransverseImpactPointExtrapolator(theField).extrapolate( impactPointState, vertexPos); //reextrapolate
192  if (!impactPointState.isValid()) return ret;
193  }
194 
195  int ndof = 2*hits.size()-5;
196  GlobalPoint vv = impactPointState.globalPosition();
197  math::XYZPoint pos( vv.x(), vv.y(), vv.z() );
198  GlobalVector pp = impactPointState.globalMomentum();
199  math::XYZVector mom( pp.x(), pp.y(), pp.z() );
200 
201  float chi2 = 0.;
202  ret = std::make_unique<reco::Track>( chi2, ndof, pos, mom,
203  impactPointState.charge(), impactPointState.curvilinearError());
204 
205 /*
206  vv = outerState.globalPosition();
207  pp = outerState.globalMomentum();
208  math::XYZPoint outerPosition( vv.x(), vv.y(), vv.z());
209  math::XYZVector outerMomentum( pp.x(), pp.y(), pp.z());
210  vv = innerState.globalPosition();
211  pp = innerState.globalMomentum();
212  math::XYZPoint innerPosition( vv.x(), vv.y(), vv.z());
213  math::XYZVector innerMomentum( pp.x(), pp.y(), pp.z());
214 
215  reco::TrackExtra extra( outerPosition, outerMomentum, true,
216  innerPosition, innerMomentum, true,
217  outerState.curvilinearError(), outerDetId,
218  innerState.curvilinearError(), innerDetId,
219  anyDirection);
220 */
221 
222 // std::cout <<"TRACK CREATED" << std::endl;
223  return ret;
224 }
float originRBound() const
bounds the particle vertex in the transverse plane
const Propagator * theOPropagator
const Propagator * thePropagator
T perp() const
Definition: PV3DBase.h:72
GlobalPoint const & origin() const
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:54
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:63
GlobalPoint globalPosition() const
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:42
const MagneticField * theField
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
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:169
const SurfaceType & surface() const
T sqr(T t)
static PlanePointer build(Args &&...args)
Definition: Plane.h:33
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
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
const edm::EventSetup * theES
const TransientTrackingRecHitBuilder * theTTRHBuilder
Definition: DetId.h:18
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
float ptMin() const
minimal pt of interest
TrajectoryStateOnSurface propagate(STA const &state, SUR const &surface) const
Definition: Propagator.h:53
TrajectoryStateOnSurface extrapolate(const FreeTrajectoryState &fts, const GlobalPoint &vtx) const
extrapolation with default (=geometrical) propagator
GlobalVector globalMomentum() const
const TrackerGeomDet * idToDet(DetId) const override
const RotationType & rotation() const
DetId geographicalId() const
T x() const
Definition: PV3DBase.h:62
const PositionType & position() const
Global3DVector GlobalVector
Definition: GlobalVector.h:10

Member Data Documentation

const reco::BeamSpot* KFBasedPixelFitter::theBeamSpot
private

Definition at line 67 of file KFBasedPixelFitter.h.

Referenced by run().

const edm::EventSetup* KFBasedPixelFitter::theES
private

Definition at line 61 of file KFBasedPixelFitter.h.

Referenced by run().

const MagneticField* KFBasedPixelFitter::theField
private

Definition at line 66 of file KFBasedPixelFitter.h.

Referenced by run().

const Propagator* KFBasedPixelFitter::theOPropagator
private

Definition at line 63 of file KFBasedPixelFitter.h.

Referenced by run().

const Propagator* KFBasedPixelFitter::thePropagator
private

Definition at line 62 of file KFBasedPixelFitter.h.

Referenced by run().

const TrackerGeometry* KFBasedPixelFitter::theTracker
private

Definition at line 65 of file KFBasedPixelFitter.h.

Referenced by run().

const TransientTrackingRecHitBuilder* KFBasedPixelFitter::theTTRHBuilder
private

Definition at line 64 of file KFBasedPixelFitter.h.

Referenced by run().