00001 #include "RecoTracker/TrackProducer/interface/GsfTrackProducerBase.h"
00002
00003 #include <memory>
00004
00005 #include "FWCore/Framework/interface/Frameworkfwd.h"
00006 #include "FWCore/Framework/interface/MakerMacros.h"
00007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00008 #include "DataFormats/GsfTrackReco/interface/GsfTrackExtra.h"
00009 #include "DataFormats/GsfTrackReco/interface/GsfTrackFwd.h"
00010 #include "TrackingTools/GeomPropagators/interface/Propagator.h"
00011 #include "TrackingTools/GeomPropagators/interface/AnalyticalPropagator.h"
00012 #include "TrackingTools/GsfTools/interface/GsfPropagatorAdapter.h"
00013 #include "TrackingTools/GsfTools/interface/MultiGaussianStateTransform.h"
00014 #include "TrackingTools/GsfTools/interface/MultiGaussianState1D.h"
00015 #include "TrackingTools/GsfTools/interface/GaussianSumUtilities1D.h"
00016 #include "TrackingTools/PatternTools/interface/TransverseImpactPointExtrapolator.h"
00017 #include "TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h"
00018 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00019 #include "TrackingTools/PatternTools/interface/CollinearFitAtTM.h"
00020
00021 #include "TrackingTools/GsfTracking/interface/TrajGsfTrackAssociation.h"
00022
00023 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00024 #include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
00025
00026 void
00027 GsfTrackProducerBase::putInEvt(edm::Event& evt,
00028 const Propagator* prop,
00029 const MeasurementTracker* measTk,
00030 std::auto_ptr<TrackingRecHitCollection>& selHits,
00031 std::auto_ptr<reco::GsfTrackCollection>& selTracks,
00032 std::auto_ptr<reco::TrackExtraCollection>& selTrackExtras,
00033 std::auto_ptr<reco::GsfTrackExtraCollection>& selGsfTrackExtras,
00034 std::auto_ptr<std::vector<Trajectory> >& selTrajectories,
00035 AlgoProductCollection& algoResults,
00036 const reco::BeamSpot& bs)
00037 {
00038
00039 TrackingRecHitRefProd rHits = evt.getRefBeforePut<TrackingRecHitCollection>();
00040 reco::TrackExtraRefProd rTrackExtras = evt.getRefBeforePut<reco::TrackExtraCollection>();
00041 reco::GsfTrackExtraRefProd rGsfTrackExtras = evt.getRefBeforePut<reco::GsfTrackExtraCollection>();
00042 reco::GsfTrackRefProd rTracks = evt.getRefBeforePut<reco::GsfTrackCollection>();
00043
00044 edm::Ref<reco::TrackExtraCollection>::key_type idx = 0;
00045 edm::Ref<reco::TrackExtraCollection>::key_type hidx = 0;
00046 edm::Ref<reco::GsfTrackExtraCollection>::key_type idxGsf = 0;
00047 edm::Ref<reco::GsfTrackCollection>::key_type iTkRef = 0;
00048 edm::Ref< std::vector<Trajectory> >::key_type iTjRef = 0;
00049 std::map<unsigned int, unsigned int> tjTkMap;
00050
00051 TSCBLBuilderNoMaterial tscblBuilder;
00052
00053 for(AlgoProductCollection::iterator i=algoResults.begin(); i!=algoResults.end();i++){
00054 Trajectory * theTraj = (*i).first;
00055 if(trajectoryInEvent_) {
00056 selTrajectories->push_back(*theTraj);
00057 iTjRef++;
00058 }
00059
00060
00061 TrajectoryFitter::RecHitContainer transHits = theTraj->recHits(useSplitting);
00062 reco::GsfTrack * theTrack = (*i).second.first;
00063 PropagationDirection seedDir = (*i).second.second;
00064
00065 LogDebug("TrackProducer") << "In GsfTrackProducerBase::putInEvt - seedDir=" << seedDir;
00066
00067 reco::GsfTrack t = * theTrack;
00068 selTracks->push_back( t );
00069 iTkRef++;
00070
00071
00072 if(trajectoryInEvent_) tjTkMap[iTjRef-1] = iTkRef-1;
00073
00074
00075 TrajectoryStateOnSurface outertsos;
00076 TrajectoryStateOnSurface innertsos;
00077 unsigned int innerId, outerId;
00078
00079
00080
00081 if (theTraj->direction() == alongMomentum) {
00082 outertsos = theTraj->lastMeasurement().updatedState();
00083 innertsos = theTraj->firstMeasurement().updatedState();
00084 outerId = theTraj->lastMeasurement().recHit()->geographicalId().rawId();
00085 innerId = theTraj->firstMeasurement().recHit()->geographicalId().rawId();
00086 } else {
00087 outertsos = theTraj->firstMeasurement().updatedState();
00088 innertsos = theTraj->lastMeasurement().updatedState();
00089 outerId = theTraj->firstMeasurement().recHit()->geographicalId().rawId();
00090 innerId = theTraj->lastMeasurement().recHit()->geographicalId().rawId();
00091 }
00092
00093 GlobalPoint v = outertsos.globalParameters().position();
00094 GlobalVector p = outertsos.globalParameters().momentum();
00095 math::XYZVector outmom( p.x(), p.y(), p.z() );
00096 math::XYZPoint outpos( v.x(), v.y(), v.z() );
00097 v = innertsos.globalParameters().position();
00098 p = innertsos.globalParameters().momentum();
00099 math::XYZVector inmom( p.x(), p.y(), p.z() );
00100 math::XYZPoint inpos( v.x(), v.y(), v.z() );
00101
00102 reco::TrackExtraRef teref= reco::TrackExtraRef ( rTrackExtras, idx ++ );
00103 reco::GsfTrack & track = selTracks->back();
00104 track.setExtra( teref );
00105
00106
00107 if (theSchool.isValid())
00108 {
00109 NavigationSetter setter( *theSchool );
00110 setSecondHitPattern(theTraj,track,prop,measTk);
00111 }
00112
00113
00114 selTrackExtras->push_back( reco::TrackExtra (outpos, outmom, true, inpos, inmom, true,
00115 outertsos.curvilinearError(), outerId,
00116 innertsos.curvilinearError(), innerId,
00117 seedDir, theTraj->seedRef()));
00118
00119
00120 reco::TrackExtra & tx = selTrackExtras->back();
00121
00122
00123 size_t ih = 0;
00124
00125
00126 if (theTraj->direction() == alongMomentum) {
00127 for( TrajectoryFitter::RecHitContainer::const_iterator j = transHits.begin();
00128 j != transHits.end(); j ++ ) {
00129 if ((**j).hit()!=0){
00130 TrackingRecHit * hit = (**j).hit()->clone();
00131 track.setHitPattern( * hit, ih ++ );
00132 selHits->push_back( hit );
00133 tx.add( TrackingRecHitRef( rHits, hidx ++ ) );
00134 }
00135 }
00136 }else{
00137 for( TrajectoryFitter::RecHitContainer::const_iterator j = transHits.end()-1;
00138 j != transHits.begin()-1; --j ) {
00139 if ((**j).hit()!=0){
00140 TrackingRecHit * hit = (**j).hit()->clone();
00141 track.setHitPattern( * hit, ih ++ );
00142 selHits->push_back( hit );
00143 tx.add( TrackingRecHitRef( rHits, hidx ++ ) );
00144 }
00145 }
00146 }
00147
00148
00149 std::vector<reco::GsfTangent> tangents;
00150 const Trajectory::DataContainer& measurements = theTraj->measurements();
00151 if ( measurements.size()>2 ) {
00152 tangents.reserve(measurements.size()-2);
00153 Trajectory::DataContainer::const_iterator ibegin,iend;
00154 int increment(0);
00155 if (theTraj->direction() == alongMomentum) {
00156 ibegin = measurements.begin() + 1;
00157 iend = measurements.end() - 1;
00158 increment = 1;
00159 }
00160 else {
00161 ibegin = measurements.end() - 2;
00162 iend = measurements.begin();
00163 increment = -1;
00164 }
00165 math::XYZPoint position;
00166 math::XYZVector momentum;
00167 Measurement1D deltaP;
00168
00169 for ( Trajectory::DataContainer::const_iterator i=ibegin;
00170 i!=iend; i+=increment ) {
00171 if ( i->recHit().get() ) {
00172 DetId detId(i->recHit()->geographicalId());
00173 if ( detId.det()==DetId::Tracker ) {
00174 int subdetId = detId.subdetId();
00175 if ( subdetId==SiStripDetId::TIB || subdetId==SiStripDetId::TID ||
00176 subdetId==SiStripDetId::TOB || subdetId==SiStripDetId::TEC ) {
00177 if ( SiStripDetId(detId).stereo() ) continue;
00178 }
00179 }
00180 }
00181 bool valid = computeModeAtTM(*i,position,momentum,deltaP);
00182 if ( valid ) {
00183 tangents.push_back(reco::GsfTangent(position,momentum,deltaP));
00184 }
00185 }
00186 }
00187
00188
00189
00190 std::vector<reco::GsfComponent5D> outerStates;
00191 outerStates.reserve(outertsos.components().size());
00192 fillStates(outertsos,outerStates);
00193 std::vector<reco::GsfComponent5D> innerStates;
00194 innerStates.reserve(innertsos.components().size());
00195 fillStates(innertsos,innerStates);
00196
00197
00198 reco::GsfTrackExtraRef terefGsf = reco::GsfTrackExtraRef ( rGsfTrackExtras, idxGsf ++ );
00199 track.setGsfExtra( terefGsf );
00200 selGsfTrackExtras->push_back( reco::GsfTrackExtra (outerStates, outertsos.localParameters().pzSign(),
00201 innerStates, innertsos.localParameters().pzSign(),
00202 tangents));
00203
00204 if ( innertsos.isValid() ) {
00205 GsfPropagatorAdapter gsfProp(AnalyticalPropagator(innertsos.magneticField(),anyDirection));
00206 TransverseImpactPointExtrapolator tipExtrapolator(gsfProp);
00207 fillMode(track,innertsos,gsfProp,tipExtrapolator,tscblBuilder,bs);
00208 }
00209
00210 delete theTrack;
00211 delete theTraj;
00212 }
00213
00214 LogTrace("TrackingRegressionTest") << "========== TrackProducer Info ===================";
00215 LogTrace("TrackingRegressionTest") << "number of finalGsfTracks: " << selTracks->size();
00216 for (reco::GsfTrackCollection::const_iterator it = selTracks->begin(); it != selTracks->end(); it++) {
00217 LogTrace("TrackingRegressionTest") << "track's n valid and invalid hit, chi2, pt : "
00218 << it->found() << " , "
00219 << it->lost() <<" , "
00220 << it->normalizedChi2() << " , "
00221 << it->pt() << " , "
00222 << it->eta() ;
00223 }
00224 LogTrace("TrackingRegressionTest") << "=================================================";
00225
00226
00227 rTracks_ = evt.put( selTracks );
00228 evt.put( selTrackExtras );
00229 evt.put( selGsfTrackExtras );
00230 evt.put( selHits );
00231
00232 if(trajectoryInEvent_) {
00233 edm::OrphanHandle<std::vector<Trajectory> > rTrajs = evt.put(selTrajectories);
00234
00235
00236 std::auto_ptr<TrajGsfTrackAssociationCollection> trajTrackMap( new TrajGsfTrackAssociationCollection() );
00237 for ( std::map<unsigned int, unsigned int>::iterator i = tjTkMap.begin();
00238 i != tjTkMap.end(); i++ ) {
00239 edm::Ref<std::vector<Trajectory> > trajRef( rTrajs, (*i).first );
00240 edm::Ref<reco::GsfTrackCollection> tkRef( rTracks_, (*i).second );
00241 trajTrackMap->insert( edm::Ref<std::vector<Trajectory> >( rTrajs, (*i).first ),
00242 edm::Ref<reco::GsfTrackCollection>( rTracks_, (*i).second ) );
00243 }
00244 evt.put( trajTrackMap );
00245 }
00246 }
00247
00248 void
00249 GsfTrackProducerBase::fillStates (TrajectoryStateOnSurface tsos,
00250 std::vector<reco::GsfComponent5D>& states) const
00251 {
00252 reco::GsfComponent5D::ParameterVector pLocS;
00253 reco::GsfComponent5D::CovarianceMatrix cLocS;
00254 std::vector<TrajectoryStateOnSurface> components(tsos.components());
00255 for ( std::vector<TrajectoryStateOnSurface>::const_iterator i=components.begin();
00256 i!=components.end(); ++i ) {
00257 states.push_back(reco::GsfComponent5D(i->weight(),i->localParameters().vector(),i->localError().matrix()));
00258 }
00259 }
00260
00261 void
00262 GsfTrackProducerBase::fillMode (reco::GsfTrack& track, const TrajectoryStateOnSurface innertsos,
00263 const Propagator& gsfProp,
00264 const TransverseImpactPointExtrapolator& tipExtrapolator,
00265 TrajectoryStateClosestToBeamLineBuilder& tscblBuilder,
00266 const reco::BeamSpot& bs) const
00267 {
00268
00269
00270
00271 GlobalPoint bsPos(bs.position().x()+(track.vz()-bs.position().z())*bs.dxdz(),
00272 bs.position().y()+(track.vz()-bs.position().z())*bs.dydz(),
00273 track.vz());
00274 TrajectoryStateOnSurface vtxTsos = tipExtrapolator.extrapolate(innertsos,bsPos);
00275 if ( !vtxTsos.isValid() ) vtxTsos = innertsos;
00276
00277 vtxTsos = gsfProp.propagate(innertsos,vtxTsos.surface());
00278 if ( !vtxTsos.isValid() ) return;
00279
00280
00281 AlgebraicVector5 modeParameters;
00282 AlgebraicSymMatrix55 modeCovariance;
00283
00284 for ( unsigned int iv=0; iv<5; ++iv ) {
00285 MultiGaussianState1D state1D = MultiGaussianStateTransform::multiState1D(vtxTsos,iv);
00286 GaussianSumUtilities1D utils(state1D);
00287 modeParameters(iv) = utils.mode().mean();
00288 modeCovariance(iv,iv) = utils.mode().variance();
00289 if ( !utils.modeIsValid() ) {
00290
00291 modeParameters(iv) = utils.mean();
00292 modeCovariance(iv,iv) = utils.variance();
00293 }
00294 }
00295
00296
00297 const AlgebraicSymMatrix55& meanCovariance(vtxTsos.localError().matrix());
00298 for ( unsigned int iv1=0; iv1<5; ++iv1 ) {
00299 for ( unsigned int iv2=0; iv2<iv1; ++iv2 ) {
00300 double cov12 = meanCovariance(iv1,iv2) *
00301 sqrt(modeCovariance(iv1,iv1)/meanCovariance(iv1,iv1)*
00302 modeCovariance(iv2,iv2)/meanCovariance(iv2,iv2));
00303 modeCovariance(iv1,iv2) = modeCovariance(iv2,iv1) = cov12;
00304 }
00305 }
00306 TrajectoryStateOnSurface modeTsos(LocalTrajectoryParameters(modeParameters,
00307 vtxTsos.localParameters().pzSign()),
00308 LocalTrajectoryError(modeCovariance),
00309 vtxTsos.surface(),
00310 vtxTsos.magneticField(),
00311 vtxTsos.surfaceSide());
00312 TrajectoryStateClosestToBeamLine tscbl = tscblBuilder(*modeTsos.freeState(),bs);
00313 if ( !tscbl.isValid() ) return;
00314
00315
00316
00317 FreeTrajectoryState fts = tscbl.trackStateAtPCA();
00318 GlobalVector tscblMom = fts.momentum();
00319 reco::GsfTrack::Vector mom(tscblMom.x(),tscblMom.y(),tscblMom.z());
00320 reco::GsfTrack::CovarianceMatrixMode cov;
00321 const AlgebraicSymMatrix55& tscblCov = fts.curvilinearError().matrix();
00322 for ( unsigned int iv1=0; iv1<reco::GsfTrack::dimensionMode; ++iv1 ) {
00323 for ( unsigned int iv2=0; iv2<reco::GsfTrack::dimensionMode; ++iv2 ) {
00324 cov(iv1,iv2) = tscblCov(iv1,iv2);
00325 }
00326 }
00327 track.setMode(fts.charge(),mom,cov);
00328 }
00329
00330 void
00331 GsfTrackProducerBase::localParametersFromQpMode (const TrajectoryStateOnSurface tsos,
00332 AlgebraicVector5& parameters,
00333 AlgebraicSymMatrix55& covariance) const
00334 {
00335
00336
00337
00338 parameters = tsos.localParameters().vector();
00339 covariance = tsos.localError().matrix();
00340
00341
00342
00343 MultiGaussianState1D qpState(MultiGaussianStateTransform::multiState1D(tsos,0));
00344 GaussianSumUtilities1D qpGS(qpState);
00345 if ( !qpGS.modeIsValid() ) return;
00346 double qp = qpGS.mode().mean();
00347 double varQp = qpGS.mode().variance();
00348
00349
00350
00351
00352 double VarQpRatio = sqrt(varQp/covariance(0,0));
00353 parameters(0) = qp;
00354 covariance(0,0) = varQp;
00355 for ( int i=1; i<5; ++i ) covariance(i,0) *= VarQpRatio;
00356 }
00357
00358 bool
00359 GsfTrackProducerBase::computeModeAtTM (const TrajectoryMeasurement& tm,
00360 reco::GsfTrackExtra::Point& position,
00361 reco::GsfTrackExtra::Vector& momentum,
00362 Measurement1D& deltaP) const
00363 {
00364
00365
00366
00367 TrajectoryStateOnSurface fwdState = tm.forwardPredictedState();
00368 TrajectoryStateOnSurface bwdState = tm.backwardPredictedState();
00369 TrajectoryStateOnSurface upState = tm.updatedState();
00370
00371
00372 if ( !fwdState.isValid() || !bwdState.isValid() || !upState.isValid() ) {
00373 return false;
00374 }
00375
00376
00377
00378
00379 GlobalPoint pos = upState.globalPosition();
00380 position = reco::GsfTrackExtra::Point(pos.x(),pos.y(),pos.z());
00381 MultiTrajectoryStateMode mts;
00382 GlobalVector mom;
00383 bool result = mts.momentumFromModeCartesian(upState,mom);
00384 if ( !result ) {
00385
00386 return false;
00387 }
00388 momentum = reco::GsfTrackExtra::Vector(mom.x(),mom.y(),mom.z());
00389
00390
00391
00392
00393
00394 AlgebraicVector5 fwdPars = fwdState.localParameters().vector();
00395 AlgebraicSymMatrix55 fwdCov = fwdState.localError().matrix();
00396 localParametersFromQpMode(fwdState,fwdPars,fwdCov);
00397 AlgebraicVector5 bwdPars = bwdState.localParameters().vector();
00398 AlgebraicSymMatrix55 bwdCov = bwdState.localError().matrix();
00399 localParametersFromQpMode(bwdState,bwdPars,bwdCov);
00400 LocalPoint hitPos(0.,0.,0.);
00401 LocalError hitErr(-1.,-1.,-1.);
00402 if ( tm.recHit()->isValid() ) {
00403 hitPos = tm.recHit()->localPosition();
00404 hitErr = tm.recHit()->localPositionError();
00405 }
00406 CollinearFitAtTM2 collinearFit(fwdPars,fwdCov,bwdPars,bwdCov,hitPos,hitErr);
00407 deltaP = collinearFit.deltaP();
00408
00409 return true;
00410 }