CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
SimpleTrackRefitter Class Reference

#include <SimpleTrackRefitter.h>

Public Member Functions

std::vector< TrajectoryrefitTrack (const reco::Track &newTrack, const uint32_t ExcludedDetId=0)
 The main methods. More...
 
std::vector< TrajectoryrefitTrack (const TrajectorySeed &seed, const TrackingRecHitCollection &hits, const uint32_t ExcludedDetid=0)
 
std::vector< TrajectoryrefitTrack (const TrajectorySeed &seed, const reco::Track &theT, const uint32_t ExcludedDetid=0)
 
void setServices (const edm::EventSetup &es)
 
 SimpleTrackRefitter (const edm::ParameterSet &)
 Constructor. More...
 
virtual ~SimpleTrackRefitter ()
 Destructor. More...
 

Private Member Functions

Trajectory createStartingTrajectory (const TrajectorySeed &seed) const
 
void destroyServices ()
 
void initServices (const bool &seedAlongMomentum)
 
std::vector
< TrajectoryMeasurement
seedMeasurements (const TrajectorySeed &seed) const
 
TrajectoryStateOnSurface startingTSOS (const TrajectorySeed &seed) const
 

Private Attributes

edm::ParameterSet conf_
 
bool isCosmics_
 
const MagneticFieldmagfield
 
const
TransientTrackingRecHitBuilder
RHBuilder
 
const Chi2MeasurementEstimatortheEstimator
 
const TrajectoryFittertheFitter
 
const PropagatorthePropagator
 
const PropagatorthePropagatorOp
 
const TrajectorySmoothertheSmoother
 
const KFUpdatortheUpdator
 
const TrackerGeometrytracker
 
const TrajectoryStateTransform tsTransform
 

Detailed Description

This class refits a track (from cosmic or standard tracking). The final result is a Trajectory refitted and smoothed. To make the refitting (and the smoothing) the usual KF tools are used. An arbitrary DetUnit can be ignored during the fit, so that the resulting trajectory can be used for unbiased residual studies on that Det.

Definition at line 46 of file SimpleTrackRefitter.h.

Constructor & Destructor Documentation

SimpleTrackRefitter::SimpleTrackRefitter ( const edm::ParameterSet iConfig)

Constructor.

Definition at line 25 of file SimpleTrackRefitter.cc.

References conf_, edm::ParameterSet::getParameter(), and isCosmics_.

25  :
27 {
28  isCosmics_ = conf_.getParameter<bool>("cosmic");
29 }
T getParameter(std::string const &) const
const TransientTrackingRecHitBuilder * RHBuilder
edm::ParameterSet conf_
const KFUpdator * theUpdator
const TrajectoryStateTransform tsTransform
const Propagator * thePropagator
const TrajectoryFitter * theFitter
const Chi2MeasurementEstimator * theEstimator
const Propagator * thePropagatorOp
const TrajectorySmoother * theSmoother
SimpleTrackRefitter::~SimpleTrackRefitter ( )
virtual

Destructor.

Definition at line 31 of file SimpleTrackRefitter.cc.

32 {
33 
34 }

Member Function Documentation

Trajectory SimpleTrackRefitter::createStartingTrajectory ( const TrajectorySeed seed) const
private

Definition at line 251 of file SimpleTrackRefitter.cc.

References TrajectorySeed::direction(), i, query::result, and seedMeasurements().

Referenced by refitTrack().

252 {
253  Trajectory result( seed, seed.direction());
254  std::vector<TrajectoryMeasurement> seedMeas = seedMeasurements(seed);
255  if ( !seedMeas.empty()) {
256  for (std::vector<TrajectoryMeasurement>::const_iterator i=seedMeas.begin(); i!=seedMeas.end(); i++){
257  result.push(*i);
258  }
259  }
260  return result;
261 }
PropagationDirection direction() const
int i
Definition: DBlmapReader.cc:9
tuple result
Definition: query.py:137
std::vector< TrajectoryMeasurement > seedMeasurements(const TrajectorySeed &seed) const
void SimpleTrackRefitter::destroyServices ( )
private

Definition at line 88 of file SimpleTrackRefitter.cc.

References theEstimator, theFitter, thePropagator, thePropagatorOp, theSmoother, and theUpdator.

Referenced by refitTrack().

89 {
90  if(thePropagator) { delete thePropagator; thePropagator=0; }
92  if(theUpdator) { delete theUpdator; theUpdator=0; }
93  if(theEstimator) { delete theEstimator; theEstimator=0; }
94  if(theSmoother) { delete theSmoother; theSmoother=0; }
95  if(theFitter) { delete theFitter; theFitter=0; }
96 }
const KFUpdator * theUpdator
const Propagator * thePropagator
const TrajectoryFitter * theFitter
const Chi2MeasurementEstimator * theEstimator
const Propagator * thePropagatorOp
const TrajectorySmoother * theSmoother
void SimpleTrackRefitter::initServices ( const bool &  seedAlongMomentum)
private

Definition at line 36 of file SimpleTrackRefitter.cc.

References alongMomentum, Chi2MeasurementEstimatorESProducer_cfi::Chi2MeasurementEstimator, KFTrajectoryFitterESProducer_cfi::KFTrajectoryFitter, KFTrajectorySmootherESProducer_cfi::KFTrajectorySmoother, magfield, oppositeToMomentum, theEstimator, theFitter, thePropagator, thePropagatorOp, theSmoother, and theUpdator.

Referenced by refitTrack().

37 {
38  // dynamic services for cosmic refit
39  if(seedAlongMomentum) {
42  }
43  else {
46  }
47  theUpdator= new KFUpdator();
50  *theUpdator,
51  *theEstimator) ;
53  *theUpdator,
54  *theEstimator);
55 }
const MagneticField * magfield
const KFUpdator * theUpdator
const Propagator * thePropagator
const TrajectoryFitter * theFitter
const Chi2MeasurementEstimator * theEstimator
const Propagator * thePropagatorOp
const TrajectorySmoother * theSmoother
std::vector< Trajectory > SimpleTrackRefitter::refitTrack ( const reco::Track newTrack,
const uint32_t  ExcludedDetId = 0 
)

The main methods.

Definition at line 98 of file SimpleTrackRefitter.cc.

References TransientTrackingRecHitBuilder::build(), funct::C, TrajectoryFitter::fit(), i, TrackingRecHit::inactive, LogDebug, Propagator::magneticField(), NULL, Propagator::propagate(), reco::Track::recHitsBegin(), reco::Track::recHitsEnd(), RHBuilder, theFitter, thePropagator, tmp, and DetId::Tracker.

Referenced by SiStripFineDelayTLA::findtrackangle().

99 {
100  // convert the TrackingRecHit vector to a TransientTrackingRecHit vector
101  LogDebug("SimpleTrackRefitter") << "Start\n";
103  for (trackingRecHit_iterator i=theT.recHitsBegin();i!=theT.recHitsEnd(); i++) {
104  if ( (*i)->geographicalId().det() == DetId::Tracker ) {
105  if((*i)->isValid()) {
106  if((*i)->geographicalId().rawId()!=ExcludedDetid) {
107  tmp.push_back(RHBuilder->build(&**i ));
108  } else {
109  InvalidTrackingRecHit* irh = new InvalidTrackingRecHit((*i)->geographicalId(), TrackingRecHit::inactive);
110  tmp.push_back(RHBuilder->build(irh));
111  delete irh;
112  }
113  } else {
114  tmp.push_back(RHBuilder->build(&**i ));
115  }
116  }
117  }
118  LogDebug("SimpleTrackRefitter") << "Transient rechit filled" << "\n";
119  // sort rechits alongmomentum
122  for (TransientTrackingRecHit::RecHitContainer::const_iterator it=tmp.begin(); it!=tmp.end();it++){
123  if ((**it).isValid()) {
124  firstHit = &(*it);
125  break;
126  }
127  }
129  for (TransientTrackingRecHit::RecHitContainer::const_iterator it=tmp.end()-1; it!=tmp.begin()-1;it--){
130  if ((**it).isValid()) {
131  lastHit= &(*it);
132  break;
133  }
134  }
135  if ((*firstHit)->globalPosition().mag2() > ((*lastHit)->globalPosition().mag2()) ){
136  for (TransientTrackingRecHit::RecHitContainer::const_iterator it=tmp.end()-1;it!=tmp.begin()-1;it--){
137  hits.push_back(*it);
138  }
139  } else hits=tmp;
140 
141  // build the transient track and fit it
142  std::vector<Trajectory> trajVec;
144  TrajectoryStateOnSurface firstState=thePropagator->propagate(theTT.impactPointState(), hits.front()->det()->surface());
145  AlgebraicSymMatrix C(5,1);
146  C *= 100.;
147  if(!firstState.isValid()) return trajVec;
148  TrajectoryStateOnSurface theTSOS( firstState.localParameters(), LocalTrajectoryError(C),
149  firstState.surface(), thePropagator->magneticField());
150  LogDebug("SimpleTrackRefitter") << "initial TSOS\n" << theTSOS << "\n";
151  const TrajectorySeed seed;//empty seed: not needed
152  trajVec = theFitter->fit(seed, hits, theTSOS);
153  LogDebug("SimpleTrackRefitter") <<" FITTER FOUND "<< trajVec.size() << " TRAJECTORIES" <<"\n";
154  return trajVec;
155 }
#define LogDebug(id)
const TransientTrackingRecHitBuilder * RHBuilder
int i
Definition: DBlmapReader.cc:9
#define NULL
Definition: scimark2.h:8
std::vector< ConstRecHitPointer > RecHitContainer
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
const Propagator * thePropagator
virtual TrajectoryStateOnSurface propagate(const FreeTrajectoryState &, const Surface &) const
Definition: Propagator.cc:9
const TrajectoryFitter * theFitter
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
virtual const MagneticField * magneticField() const =0
virtual std::vector< Trajectory > fit(const Trajectory &) const =0
CLHEP::HepSymMatrix AlgebraicSymMatrix
std::vector< Trajectory > SimpleTrackRefitter::refitTrack ( const TrajectorySeed seed,
const TrackingRecHitCollection hits,
const uint32_t  ExcludedDetid = 0 
)

Definition at line 157 of file SimpleTrackRefitter.cc.

References alongMomentum, TransientTrackingRecHitBuilder::build(), createStartingTrajectory(), destroyServices(), TrajectorySeed::direction(), Chi2MeasurementEstimator::estimate(), TrajectoryFitter::fit(), TrackingRecHit::geographicalId(), TrackerGeometry::idToDet(), TrackingRecHit::inactive, initServices(), TrajectoryStateOnSurface::isValid(), Propagator::propagate(), DetId::rawId(), RHBuilder, edm::OwnVector< T, P >::size(), GeomDet::surface(), theEstimator, theFitter, thePropagator, thePropagatorOp, theSmoother, theUpdator, tracker, TrajectorySmoother::trajectories(), and KFUpdator::update().

160 {
164  for (unsigned int icosmhit=hits.size()-1;icosmhit+1>0;icosmhit--){
166  const TrackingRecHit* rh = &(hits[icosmhit]);
167  if(rh->geographicalId().rawId() == ExcludedDetid) {
169  tmphit=RHBuilder->build(irh);
170  delete irh;
171  }
172  else {
173  tmphit=RHBuilder->build(rh);
174  }
175  trans_hits.push_back(&(*tmphit));
176  if (icosmhit<hits.size()-1){
177  TSOS prSt = thePropagator->propagate(traj.lastMeasurement().updatedState(),
178  tracker->idToDet(hits[icosmhit].geographicalId())->surface());
179  if(prSt.isValid() && tmphit->isValid()) {
180  if(theUpdator->update(prSt,*tmphit).isValid()) {
181  traj.push(TrajectoryMeasurement(prSt, theUpdator->update(prSt,*tmphit),
182  tmphit, theEstimator->estimate(prSt,*tmphit).second));
183  }
184  }
185  }
186  }
187  std::vector<Trajectory> smoothtraj;
188  if (thePropagatorOp->propagate(traj.lastMeasurement().updatedState(),
189  tracker->idToDet((*trans_hits.begin())->geographicalId())->surface()).isValid()){
190  TSOS startingState= TrajectoryStateWithArbitraryError()
191  (thePropagatorOp->propagate(traj.lastMeasurement().updatedState(),
192  tracker->idToDet((*trans_hits.begin())->geographicalId())->surface()));
193  std::vector<Trajectory> fittraj=theFitter->fit(seed,trans_hits,startingState);
194  if (fittraj.size()) smoothtraj=theSmoother->trajectories(*(fittraj.begin()));
195  }
196  destroyServices();
197  return smoothtraj;
198 }
PropagationDirection direction() const
const TransientTrackingRecHitBuilder * RHBuilder
void initServices(const bool &seedAlongMomentum)
size_type size() const
Definition: OwnVector.h:262
std::vector< ConstRecHitPointer > RecHitContainer
virtual TrajectoryContainer trajectories(const Trajectory &) const =0
const TrackerGeometry * tracker
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
const KFUpdator * theUpdator
const Propagator * thePropagator
virtual const GeomDet * idToDet(DetId) const
virtual std::pair< bool, double > estimate(const TrajectoryStateOnSurface &, const TransientTrackingRecHit &) const
Trajectory createStartingTrajectory(const TrajectorySeed &seed) const
virtual TrajectoryStateOnSurface propagate(const FreeTrajectoryState &, const Surface &) const
Definition: Propagator.cc:9
const TrajectoryFitter * theFitter
bool isValid() const
const Chi2MeasurementEstimator * theEstimator
virtual std::vector< Trajectory > fit(const Trajectory &) const =0
const Propagator * thePropagatorOp
TrajectoryStateOnSurface update(const TrajectoryStateOnSurface &, const TransientTrackingRecHit &) const
Definition: KFUpdator.cc:10
DetId geographicalId() const
virtual const BoundPlane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
const TrajectorySmoother * theSmoother
std::vector< Trajectory > SimpleTrackRefitter::refitTrack ( const TrajectorySeed seed,
const reco::Track theT,
const uint32_t  ExcludedDetid = 0 
)

Definition at line 200 of file SimpleTrackRefitter.cc.

References alongMomentum, TransientTrackingRecHitBuilder::build(), createStartingTrajectory(), destroyServices(), TrajectorySeed::direction(), Chi2MeasurementEstimator::estimate(), TrajectoryFitter::fit(), i, TrackerGeometry::idToDet(), TrackingRecHit::inactive, initServices(), TrajectoryStateOnSurface::isValid(), Propagator::propagate(), reco::Track::recHitsBegin(), reco::Track::recHitsEnd(), RHBuilder, theEstimator, theFitter, thePropagator, thePropagatorOp, theSmoother, theUpdator, tracker, TrajectorySmoother::trajectories(), and KFUpdator::update().

203 {
207  for (trackingRecHit_iterator i=theT.recHitsBegin();i!=theT.recHitsEnd(); i++) {
209  if((*i)->geographicalId().rawId() == ExcludedDetid) {
210  InvalidTrackingRecHit* irh = new InvalidTrackingRecHit((*i)->geographicalId(), TrackingRecHit::inactive);
211  tmphit=RHBuilder->build(irh);
212  delete irh;
213  }
214  else {
215  tmphit=RHBuilder->build(&**i);
216  }
217  trans_hits.push_back(&(*tmphit));
218  if (i!=theT.recHitsBegin()){
219  TSOS prSt = thePropagator->propagate(traj.lastMeasurement().updatedState(),
220  tracker->idToDet((*i)->geographicalId())->surface());
221  if(prSt.isValid() && tmphit->isValid()) {
222  if(theUpdator->update(prSt,*tmphit).isValid()) {
223  traj.push(TrajectoryMeasurement(prSt, theUpdator->update(prSt,*tmphit),
224  tmphit, theEstimator->estimate(prSt,*tmphit).second));
225  }
226  }
227  }
228  }
229  std::vector<Trajectory> smoothtraj;
230  if (thePropagatorOp->propagate(traj.lastMeasurement().updatedState(),
231  tracker->idToDet((*trans_hits.begin())->geographicalId())->surface()).isValid()){
232  TSOS startingState= TrajectoryStateWithArbitraryError()
233  (thePropagatorOp->propagate(traj.lastMeasurement().updatedState(),
234  tracker->idToDet((*trans_hits.begin())->geographicalId())->surface()));
235  std::vector<Trajectory> fittraj=theFitter->fit(seed,trans_hits,startingState);
236  if (fittraj.size()) smoothtraj=theSmoother->trajectories(*(fittraj.begin()));
237  }
238  destroyServices();
239  return smoothtraj;
240 }
PropagationDirection direction() const
const TransientTrackingRecHitBuilder * RHBuilder
int i
Definition: DBlmapReader.cc:9
void initServices(const bool &seedAlongMomentum)
std::vector< ConstRecHitPointer > RecHitContainer
virtual TrajectoryContainer trajectories(const Trajectory &) const =0
const TrackerGeometry * tracker
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
const KFUpdator * theUpdator
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:63
const Propagator * thePropagator
virtual const GeomDet * idToDet(DetId) const
virtual std::pair< bool, double > estimate(const TrajectoryStateOnSurface &, const TransientTrackingRecHit &) const
Trajectory createStartingTrajectory(const TrajectorySeed &seed) const
virtual TrajectoryStateOnSurface propagate(const FreeTrajectoryState &, const Surface &) const
Definition: Propagator.cc:9
const TrajectoryFitter * theFitter
bool isValid() const
const Chi2MeasurementEstimator * theEstimator
virtual std::vector< Trajectory > fit(const Trajectory &) const =0
const Propagator * thePropagatorOp
TrajectoryStateOnSurface update(const TrajectoryStateOnSurface &, const TransientTrackingRecHit &) const
Definition: KFUpdator.cc:10
const TrajectorySmoother * theSmoother
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:65
std::vector< TrajectoryMeasurement > SimpleTrackRefitter::seedMeasurements ( const TrajectorySeed seed) const
private

Definition at line 263 of file SimpleTrackRefitter.cc.

References TransientTrackingRecHitBuilder::build(), TrackerGeometry::idToDet(), TrajectorySeed::recHits(), query::result, RHBuilder, startingTSOS(), GeomDet::surface(), and tracker.

Referenced by createStartingTrajectory().

264 {
265  std::vector<TrajectoryMeasurement> result;
266  TrajectorySeed::range hitRange = seed.recHits();
267  for (TrajectorySeed::const_iterator ihit = hitRange.first; ihit != hitRange.second; ihit++) {
269  const GeomDet* hitGeomDet = tracker->idToDet( ihit->geographicalId());
270  TSOS invalidState( new BasicSingleTrajectoryState( hitGeomDet->surface()));
271  if (ihit == hitRange.second - 1) {
272  TSOS updatedState=startingTSOS(seed);
273  result.push_back(TrajectoryMeasurement( invalidState, updatedState, recHit));
274  }
275  else {
276  result.push_back(TrajectoryMeasurement( invalidState, recHit));
277  }
278  }
279  return result;
280 }
const TransientTrackingRecHitBuilder * RHBuilder
TrajectoryStateOnSurface startingTSOS(const TrajectorySeed &seed) const
const TrackerGeometry * tracker
recHitContainer::const_iterator const_iterator
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
tuple result
Definition: query.py:137
std::pair< const_iterator, const_iterator > range
virtual const GeomDet * idToDet(DetId) const
range recHits() const
virtual const BoundPlane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
void SimpleTrackRefitter::setServices ( const edm::EventSetup es)

Definition at line 57 of file SimpleTrackRefitter.cc.

References conf_, edm::EventSetup::get(), edm::ParameterSet::getParameter(), isCosmics_, LogDebug, magfield, LargeD0_PixelPairStep_cff::propagator, RHBuilder, theFitter, thePropagator, and tracker.

Referenced by SiStripFineDelayTLA::init().

58 {
59  // get geometry
61  es.get<TrackerDigiGeometryRecord>().get(estracker);
62  tracker=&(*estracker);
63  // get magnetic field
65  es.get<IdealMagneticFieldRecord>().get(esmagfield);
66  magfield=&(*esmagfield);
67  // get the builder
69  LogDebug("SimpleTrackRefitter") << "get also the TransientTrackingRecHitBuilder" << "\n";
70  std::string builderName = conf_.getParameter<std::string>("TTRHBuilder");
71  es.get<TransientRecHitRecord>().get(builderName,builder);
72  RHBuilder=&(*builder);
73  if(isCosmics_) return;
74  // get the fitter
76  LogDebug("SimpleTrackRefitter") << "get the fitter from the ES" << "\n";
77  std::string fitterName = conf_.getParameter<std::string>("Fitter");
78  es.get<TrajectoryFitter::Record>().get(fitterName,fitter);
79  theFitter=&(*fitter);
80  // get the propagator
82  LogDebug("SimpleTrackRefitter") << "get also the propagator" << "\n";
83  std::string propagatorName = conf_.getParameter<std::string>("Propagator");
84  es.get<TrackingComponentsRecord>().get(propagatorName,propagator);
85  thePropagator=&(*propagator);
86 }
#define LogDebug(id)
T getParameter(std::string const &) const
const TransientTrackingRecHitBuilder * RHBuilder
edm::ParameterSet conf_
const TrackerGeometry * tracker
const MagneticField * magfield
const Propagator * thePropagator
const TrajectoryFitter * theFitter
const T & get() const
Definition: EventSetup.h:55
TrajectoryStateOnSurface SimpleTrackRefitter::startingTSOS ( const TrajectorySeed seed) const
private

Definition at line 242 of file SimpleTrackRefitter.cc.

References TrackerGeometry::idToDet(), magfield, TrajectorySeed::startingState(), tracker, TrajectoryStateTransform::transientState(), and tsTransform.

Referenced by seedMeasurements().

243 {
244  PTrajectoryStateOnDet pState( seed.startingState());
245  const GeomDet* gdet = tracker->idToDet(DetId(pState.detId()));
246  TSOS State= tsTransform.transientState( pState, &(gdet->surface()), &(*magfield));
247  return State;
248 
249 }
const TrackerGeometry * tracker
const MagneticField * magfield
TrajectoryStateOnSurface transientState(const PTrajectoryStateOnDet &ts, const Surface *surface, const MagneticField *field) const
const TrajectoryStateTransform tsTransform
virtual const GeomDet * idToDet(DetId) const
Definition: DetId.h:20
PTrajectoryStateOnDet const & startingState() const

Member Data Documentation

edm::ParameterSet SimpleTrackRefitter::conf_
private

Definition at line 86 of file SimpleTrackRefitter.h.

Referenced by setServices(), and SimpleTrackRefitter().

bool SimpleTrackRefitter::isCosmics_
private

Definition at line 87 of file SimpleTrackRefitter.h.

Referenced by setServices(), and SimpleTrackRefitter().

const MagneticField* SimpleTrackRefitter::magfield
private

Definition at line 85 of file SimpleTrackRefitter.h.

Referenced by initServices(), setServices(), and startingTSOS().

const TransientTrackingRecHitBuilder* SimpleTrackRefitter::RHBuilder
private

Definition at line 80 of file SimpleTrackRefitter.h.

Referenced by refitTrack(), seedMeasurements(), and setServices().

const Chi2MeasurementEstimator* SimpleTrackRefitter::theEstimator
private

Definition at line 79 of file SimpleTrackRefitter.h.

Referenced by destroyServices(), initServices(), and refitTrack().

const TrajectoryFitter* SimpleTrackRefitter::theFitter
private

Definition at line 82 of file SimpleTrackRefitter.h.

Referenced by destroyServices(), initServices(), refitTrack(), and setServices().

const Propagator* SimpleTrackRefitter::thePropagator
private

Definition at line 76 of file SimpleTrackRefitter.h.

Referenced by destroyServices(), initServices(), refitTrack(), and setServices().

const Propagator* SimpleTrackRefitter::thePropagatorOp
private

Definition at line 77 of file SimpleTrackRefitter.h.

Referenced by destroyServices(), initServices(), and refitTrack().

const TrajectorySmoother* SimpleTrackRefitter::theSmoother
private

Definition at line 81 of file SimpleTrackRefitter.h.

Referenced by destroyServices(), initServices(), and refitTrack().

const KFUpdator* SimpleTrackRefitter::theUpdator
private

Definition at line 78 of file SimpleTrackRefitter.h.

Referenced by destroyServices(), initServices(), and refitTrack().

const TrackerGeometry* SimpleTrackRefitter::tracker
private

Definition at line 84 of file SimpleTrackRefitter.h.

Referenced by refitTrack(), seedMeasurements(), setServices(), and startingTSOS().

const TrajectoryStateTransform SimpleTrackRefitter::tsTransform
private

Definition at line 83 of file SimpleTrackRefitter.h.

Referenced by startingTSOS().