#include <DQM/SiStripCommissioningSources/plugins/tracking/SimpleTrackRefitter.h>
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.
SimpleTrackRefitter::SimpleTrackRefitter | ( | const edm::ParameterSet & | iConfig | ) |
Constructor.
Definition at line 25 of file SimpleTrackRefitter.cc.
References conf_, edm::ParameterSet::getParameter(), and isCosmics_.
00025 : 00026 thePropagator(0),thePropagatorOp(0),theUpdator(0),theEstimator(0),RHBuilder(0),theSmoother(0),theFitter(0),tsTransform(),conf_(iConfig) 00027 { 00028 isCosmics_ = conf_.getParameter<bool>("cosmic"); 00029 }
SimpleTrackRefitter::~SimpleTrackRefitter | ( | ) | [virtual] |
Trajectory SimpleTrackRefitter::createStartingTrajectory | ( | const TrajectorySeed & | seed | ) | const [private] |
Definition at line 251 of file SimpleTrackRefitter.cc.
References TrajectorySeed::direction(), i, HLT_VtxMuL3::result, and seedMeasurements().
Referenced by refitTrack().
00252 { 00253 Trajectory result( seed, seed.direction()); 00254 std::vector<TrajectoryMeasurement> seedMeas = seedMeasurements(seed); 00255 if ( !seedMeas.empty()) { 00256 for (std::vector<TrajectoryMeasurement>::const_iterator i=seedMeas.begin(); i!=seedMeas.end(); i++){ 00257 result.push(*i); 00258 } 00259 } 00260 return result; 00261 }
void SimpleTrackRefitter::destroyServices | ( | ) | [private] |
Definition at line 88 of file SimpleTrackRefitter.cc.
References theEstimator, theFitter, thePropagator, thePropagatorOp, theSmoother, and theUpdator.
Referenced by refitTrack().
00089 { 00090 if(thePropagator) { delete thePropagator; thePropagator=0; } 00091 if(thePropagatorOp) { delete thePropagatorOp; thePropagatorOp =0; } 00092 if(theUpdator) { delete theUpdator; theUpdator=0; } 00093 if(theEstimator) { delete theEstimator; theEstimator=0; } 00094 if(theSmoother) { delete theSmoother; theSmoother=0; } 00095 if(theFitter) { delete theFitter; theFitter=0; } 00096 }
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().
00037 { 00038 // dynamic services for cosmic refit 00039 if(seedAlongMomentum) { 00040 thePropagator = new PropagatorWithMaterial(alongMomentum,0.1057,magfield ); 00041 thePropagatorOp = new PropagatorWithMaterial(oppositeToMomentum,0.1057,magfield ); 00042 } 00043 else { 00044 thePropagator = new PropagatorWithMaterial(oppositeToMomentum,0.1057,magfield ); 00045 thePropagatorOp = new PropagatorWithMaterial(alongMomentum,0.1057,magfield ); 00046 } 00047 theUpdator= new KFUpdator(); 00048 theEstimator= new Chi2MeasurementEstimator(30); 00049 theFitter= new KFTrajectoryFitter(*thePropagator, 00050 *theUpdator, 00051 *theEstimator) ; 00052 theSmoother= new KFTrajectorySmoother(*thePropagatorOp, 00053 *theUpdator, 00054 *theEstimator); 00055 }
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().
00203 { 00204 initServices(seed.direction()==alongMomentum); 00205 Trajectory traj=createStartingTrajectory(seed); 00206 TransientTrackingRecHit::RecHitContainer trans_hits; 00207 for (trackingRecHit_iterator i=theT.recHitsBegin();i!=theT.recHitsEnd(); i++) { 00208 TransientTrackingRecHit::RecHitPointer tmphit; 00209 if((*i)->geographicalId().rawId() == ExcludedDetid) { 00210 InvalidTrackingRecHit* irh = new InvalidTrackingRecHit((*i)->geographicalId(), TrackingRecHit::inactive); 00211 tmphit=RHBuilder->build(irh); 00212 delete irh; 00213 } 00214 else { 00215 tmphit=RHBuilder->build(&**i); 00216 } 00217 trans_hits.push_back(&(*tmphit)); 00218 if (i!=theT.recHitsBegin()){ 00219 TSOS prSt = thePropagator->propagate(traj.lastMeasurement().updatedState(), 00220 tracker->idToDet((*i)->geographicalId())->surface()); 00221 if(prSt.isValid() && tmphit->isValid()) { 00222 if(theUpdator->update(prSt,*tmphit).isValid()) { 00223 traj.push(TrajectoryMeasurement(prSt, theUpdator->update(prSt,*tmphit), 00224 tmphit, theEstimator->estimate(prSt,*tmphit).second)); 00225 } 00226 } 00227 } 00228 } 00229 std::vector<Trajectory> smoothtraj; 00230 if (thePropagatorOp->propagate(traj.lastMeasurement().updatedState(), 00231 tracker->idToDet((*trans_hits.begin())->geographicalId())->surface()).isValid()){ 00232 TSOS startingState= TrajectoryStateWithArbitraryError() 00233 (thePropagatorOp->propagate(traj.lastMeasurement().updatedState(), 00234 tracker->idToDet((*trans_hits.begin())->geographicalId())->surface())); 00235 std::vector<Trajectory> fittraj=theFitter->fit(seed,trans_hits,startingState); 00236 if (fittraj.size()) smoothtraj=theSmoother->trajectories(*(fittraj.begin())); 00237 } 00238 destroyServices(); 00239 return smoothtraj; 00240 }
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().
00160 { 00161 initServices(seed.direction()==alongMomentum); 00162 Trajectory traj=createStartingTrajectory(seed); 00163 TransientTrackingRecHit::RecHitContainer trans_hits; 00164 for (unsigned int icosmhit=hits.size()-1;icosmhit+1>0;icosmhit--){ 00165 TransientTrackingRecHit::RecHitPointer tmphit; 00166 const TrackingRecHit* rh = &(hits[icosmhit]); 00167 if(rh->geographicalId().rawId() == ExcludedDetid) { 00168 InvalidTrackingRecHit* irh = new InvalidTrackingRecHit(rh->geographicalId(), TrackingRecHit::inactive); 00169 tmphit=RHBuilder->build(irh); 00170 delete irh; 00171 } 00172 else { 00173 tmphit=RHBuilder->build(rh); 00174 } 00175 trans_hits.push_back(&(*tmphit)); 00176 if (icosmhit<hits.size()-1){ 00177 TSOS prSt = thePropagator->propagate(traj.lastMeasurement().updatedState(), 00178 tracker->idToDet(hits[icosmhit].geographicalId())->surface()); 00179 if(prSt.isValid() && tmphit->isValid()) { 00180 if(theUpdator->update(prSt,*tmphit).isValid()) { 00181 traj.push(TrajectoryMeasurement(prSt, theUpdator->update(prSt,*tmphit), 00182 tmphit, theEstimator->estimate(prSt,*tmphit).second)); 00183 } 00184 } 00185 } 00186 } 00187 std::vector<Trajectory> smoothtraj; 00188 if (thePropagatorOp->propagate(traj.lastMeasurement().updatedState(), 00189 tracker->idToDet((*trans_hits.begin())->geographicalId())->surface()).isValid()){ 00190 TSOS startingState= TrajectoryStateWithArbitraryError() 00191 (thePropagatorOp->propagate(traj.lastMeasurement().updatedState(), 00192 tracker->idToDet((*trans_hits.begin())->geographicalId())->surface())); 00193 std::vector<Trajectory> fittraj=theFitter->fit(seed,trans_hits,startingState); 00194 if (fittraj.size()) smoothtraj=theSmoother->trajectories(*(fittraj.begin())); 00195 } 00196 destroyServices(); 00197 return smoothtraj; 00198 }
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, it, LogDebug, Propagator::magneticField(), NULL, Propagator::propagate(), reco::Track::recHitsBegin(), reco::Track::recHitsEnd(), RHBuilder, theFitter, thePropagator, tmp, and DetId::Tracker.
Referenced by SiStripFineDelayTLA::findtrackangle().
00099 { 00100 // convert the TrackingRecHit vector to a TransientTrackingRecHit vector 00101 LogDebug("SimpleTrackRefitter") << "Start\n"; 00102 TransientTrackingRecHit::RecHitContainer tmp; 00103 for (trackingRecHit_iterator i=theT.recHitsBegin();i!=theT.recHitsEnd(); i++) { 00104 if ( (*i)->geographicalId().det() == DetId::Tracker ) { 00105 if((*i)->isValid()) { 00106 if((*i)->geographicalId().rawId()!=ExcludedDetid) { 00107 tmp.push_back(RHBuilder->build(&**i )); 00108 } else { 00109 InvalidTrackingRecHit* irh = new InvalidTrackingRecHit((*i)->geographicalId(), TrackingRecHit::inactive); 00110 tmp.push_back(RHBuilder->build(irh)); 00111 delete irh; 00112 } 00113 } else { 00114 tmp.push_back(RHBuilder->build(&**i )); 00115 } 00116 } 00117 } 00118 LogDebug("SimpleTrackRefitter") << "Transient rechit filled" << "\n"; 00119 // sort rechits alongmomentum 00120 TransientTrackingRecHit::RecHitContainer hits; 00121 const TransientTrackingRecHit::ConstRecHitPointer *firstHit = NULL; 00122 for (TransientTrackingRecHit::RecHitContainer::const_iterator it=tmp.begin(); it!=tmp.end();it++){ 00123 if ((**it).isValid()) { 00124 firstHit = &(*it); 00125 break; 00126 } 00127 } 00128 const TransientTrackingRecHit::ConstRecHitPointer *lastHit = NULL; 00129 for (TransientTrackingRecHit::RecHitContainer::const_iterator it=tmp.end()-1; it!=tmp.begin()-1;it--){ 00130 if ((**it).isValid()) { 00131 lastHit= &(*it); 00132 break; 00133 } 00134 } 00135 if ((*firstHit)->globalPosition().mag2() > ((*lastHit)->globalPosition().mag2()) ){ 00136 for (TransientTrackingRecHit::RecHitContainer::const_iterator it=tmp.end()-1;it!=tmp.begin()-1;it--){ 00137 hits.push_back(*it); 00138 } 00139 } else hits=tmp; 00140 00141 // build the transient track and fit it 00142 std::vector<Trajectory> trajVec; 00143 reco::TransientTrack theTT(theT, thePropagator->magneticField() ); 00144 TrajectoryStateOnSurface firstState=thePropagator->propagate(theTT.impactPointState(), hits.front()->det()->surface()); 00145 AlgebraicSymMatrix C(5,1); 00146 C *= 100.; 00147 if(!firstState.isValid()) return trajVec; 00148 TrajectoryStateOnSurface theTSOS( firstState.localParameters(), LocalTrajectoryError(C), 00149 firstState.surface(), thePropagator->magneticField()); 00150 LogDebug("SimpleTrackRefitter") << "initial TSOS\n" << theTSOS << "\n"; 00151 const TrajectorySeed seed;//empty seed: not needed 00152 trajVec = theFitter->fit(seed, hits, theTSOS); 00153 LogDebug("SimpleTrackRefitter") <<" FITTER FOUND "<< trajVec.size() << " TRAJECTORIES" <<"\n"; 00154 return trajVec; 00155 }
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(), HLT_VtxMuL3::result, RHBuilder, startingTSOS(), GeomDet::surface(), and tracker.
Referenced by createStartingTrajectory().
00264 { 00265 std::vector<TrajectoryMeasurement> result; 00266 TrajectorySeed::range hitRange = seed.recHits(); 00267 for (TrajectorySeed::const_iterator ihit = hitRange.first; ihit != hitRange.second; ihit++) { 00268 TransientTrackingRecHit::RecHitPointer recHit = RHBuilder->build(&(*ihit)); 00269 const GeomDet* hitGeomDet = tracker->idToDet( ihit->geographicalId()); 00270 TSOS invalidState( new BasicSingleTrajectoryState( hitGeomDet->surface())); 00271 if (ihit == hitRange.second - 1) { 00272 TSOS updatedState=startingTSOS(seed); 00273 result.push_back(TrajectoryMeasurement( invalidState, updatedState, recHit)); 00274 } 00275 else { 00276 result.push_back(TrajectoryMeasurement( invalidState, recHit)); 00277 } 00278 } 00279 return result; 00280 }
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, RHBuilder, theFitter, thePropagator, and tracker.
Referenced by SiStripFineDelayTLA::init().
00058 { 00059 // get geometry 00060 edm::ESHandle<TrackerGeometry> estracker; 00061 es.get<TrackerDigiGeometryRecord>().get(estracker); 00062 tracker=&(*estracker); 00063 // get magnetic field 00064 edm::ESHandle<MagneticField> esmagfield; 00065 es.get<IdealMagneticFieldRecord>().get(esmagfield); 00066 magfield=&(*esmagfield); 00067 // get the builder 00068 edm::ESHandle<TransientTrackingRecHitBuilder> builder; 00069 LogDebug("SimpleTrackRefitter") << "get also the TransientTrackingRecHitBuilder" << "\n"; 00070 std::string builderName = conf_.getParameter<std::string>("TTRHBuilder"); 00071 es.get<TransientRecHitRecord>().get(builderName,builder); 00072 RHBuilder=&(*builder); 00073 if(isCosmics_) return; 00074 // get the fitter 00075 edm::ESHandle<TrajectoryFitter> fitter; 00076 LogDebug("SimpleTrackRefitter") << "get the fitter from the ES" << "\n"; 00077 std::string fitterName = conf_.getParameter<std::string>("Fitter"); 00078 es.get<TrackingComponentsRecord>().get(fitterName,fitter); 00079 theFitter=&(*fitter); 00080 // get the propagator 00081 edm::ESHandle<Propagator> propagator; 00082 LogDebug("SimpleTrackRefitter") << "get also the propagator" << "\n"; 00083 std::string propagatorName = conf_.getParameter<std::string>("Propagator"); 00084 es.get<TrackingComponentsRecord>().get(propagatorName,propagator); 00085 thePropagator=&(*propagator); 00086 }
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().
00243 { 00244 PTrajectoryStateOnDet pState( seed.startingState()); 00245 const GeomDet* gdet = tracker->idToDet(DetId(pState.detId())); 00246 TSOS State= tsTransform.transientState( pState, &(gdet->surface()), &(*magfield)); 00247 return State; 00248 00249 }
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] |