00001 #include <memory>
00002
00003 #include "FWCore/Framework/interface/Event.h"
00004 #include "FWCore/Framework/interface/EventSetup.h"
00005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00006 #include "FWCore/Framework/interface/ESHandle.h"
00007
00008 #include "DataFormats/Common/interface/Handle.h"
00009 #include "DataFormats/Common/interface/OwnVector.h"
00010 #include "DataFormats/TrackCandidate/interface/TrackCandidateCollection.h"
00011 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
00012 #include "DataFormats/TrackerRecHit2D/interface/SiTrackerGSRecHit2DCollection.h"
00013 #include "DataFormats/TrackerRecHit2D/interface/SiTrackerGSMatchedRecHit2DCollection.h"
00014 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00015 #include "DataFormats/TrackReco/interface/Track.h"
00016 #include "DataFormats/TrackReco/interface/TrackExtraFwd.h"
00017
00018 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00019 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00020
00021 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00022 #include "TrackingTools/PatternTools/interface/TrajTrackAssociation.h"
00023
00024 #include "FastSimulation/Tracking/interface/TrackerRecHit.h"
00025
00026
00027 #include "FastSimulation/Tracking/plugins/TrackCandidateProducer.h"
00028
00029 #include <vector>
00030 #include <map>
00031
00032 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00033 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00034 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00035 #include "MagneticField/Engine/interface/MagneticField.h"
00036
00037 #include "SimDataFormats/Track/interface/SimTrack.h"
00038 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00039 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
00040 #include "TrackingTools/GeomPropagators/interface/Propagator.h"
00041 #include "TrackingTools/Records/interface/TrackingComponentsRecord.h"
00042
00043
00044 #include "TrackingTools/MaterialEffects/interface/PropagatorWithMaterial.h"
00045
00046 #include "SimDataFormats/Track/interface/SimTrackContainer.h"
00047 #include "SimDataFormats/Vertex/interface/SimVertexContainer.h"
00048
00049
00050
00051
00052
00053
00054 TrackCandidateProducer::TrackCandidateProducer(const edm::ParameterSet& conf):thePropagator(0)
00055 {
00056 #ifdef FAMOS_DEBUG
00057 std::cout << "TrackCandidateProducer created" << std::endl;
00058 #endif
00059
00060
00061 produces<TrackCandidateCollection>();
00062
00063
00064
00065 produces<reco::TrackCollection>();
00066 produces<TrackingRecHitCollection>();
00067 produces<reco::TrackExtraCollection>();
00068 produces<std::vector<Trajectory> >();
00069 produces<TrajTrackAssociationCollection>();
00070
00071
00072 seedProducer = conf.getParameter<edm::InputTag>("SeedProducer");
00073
00074
00075 hitProducer = conf.getParameter<edm::InputTag>("HitProducer");
00076
00077
00078
00079 trackProducers = conf.getParameter<std::vector<edm::InputTag> >("TrackProducers");
00080
00081
00082 keepFittedTracks = conf.getParameter<bool>("KeepFittedTracks");
00083
00084
00085 minNumberOfCrossedLayers = conf.getParameter<unsigned int>("MinNumberOfCrossedLayers");
00086
00087
00088 maxNumberOfCrossedLayers = conf.getParameter<unsigned int>("MaxNumberOfCrossedLayers");
00089
00090
00091 rejectOverlaps = conf.getParameter<bool>("OverlapCleaning");
00092
00093
00094 splitHits = conf.getParameter<bool>("SplitHits");
00095
00096
00097
00098 seedCleaning = conf.getParameter<bool>("SeedCleaning");
00099
00100 simTracks_ = conf.getParameter<edm::InputTag>("SimTracks");
00101 estimatorCut_= conf.getParameter<double>("EstimatorCut");
00102 }
00103
00104
00105
00106 TrackCandidateProducer::~TrackCandidateProducer() {
00107
00108 if(thePropagator) delete thePropagator;
00109
00110
00111 #ifdef FAMOS_DEBUG
00112 std::cout << "TrackCandidateProducer destructed" << std::endl;
00113 #endif
00114
00115
00116 }
00117
00118 void
00119 TrackCandidateProducer::beginRun(edm::Run const&, const edm::EventSetup & es) {
00120
00121
00122
00123 edm::ESHandle<MagneticField> magField;
00124 edm::ESHandle<TrackerGeometry> geometry;
00125
00126 es.get<IdealMagneticFieldRecord>().get(magField);
00127 es.get<TrackerDigiGeometryRecord>().get(geometry);
00128
00129 theMagField = &(*magField);
00130 theGeometry = &(*geometry);
00131
00132 thePropagator = new PropagatorWithMaterial(alongMomentum,0.105,&(*theMagField));
00133 }
00134
00135
00136 void
00137 TrackCandidateProducer::produce(edm::Event& e, const edm::EventSetup& es) {
00138
00139 #ifdef FAMOS_DEBUG
00140 std::cout << "################################################################" << std::endl;
00141 std::cout << " TrackCandidateProducer produce init " << std::endl;
00142 #endif
00143
00144
00145 typedef std::pair<reco::TrackRef,edm::Ref<std::vector<Trajectory> > > TrackPair;
00146 typedef std::map<unsigned,TrackPair> TrackMap;
00147
00148
00149 std::auto_ptr<TrackCandidateCollection> output(new TrackCandidateCollection);
00150 std::auto_ptr<reco::TrackCollection> recoTracks(new reco::TrackCollection);
00151 std::auto_ptr<TrackingRecHitCollection> recoHits(new TrackingRecHitCollection);
00152 std::auto_ptr<reco::TrackExtraCollection> recoTrackExtras(new reco::TrackExtraCollection);
00153 std::auto_ptr<std::vector<Trajectory> > recoTrajectories(new std::vector<Trajectory> );
00154 std::auto_ptr<TrajTrackAssociationCollection> recoTrajTrackMap( new TrajTrackAssociationCollection() );
00155
00156
00157
00158 edm::Handle<edm::View<TrajectorySeed> > theSeeds;
00159 e.getByLabel(seedProducer,theSeeds);
00160
00161
00162 edm::ESHandle<TrackerTopology> tTopoHand;
00163 es.get<IdealGeometryRecord>().get(tTopoHand);
00164 const TrackerTopology *tTopo=tTopoHand.product();
00165
00166
00167
00168 if(theSeeds->size() == 0) {
00169 e.put(output);
00170 e.put(recoTracks);
00171 e.put(recoHits);
00172 e.put(recoTrackExtras);
00173 e.put(recoTrajectories);
00174 e.put(recoTrajTrackMap);
00175 return;
00176 }
00177
00178
00179
00180 edm::Handle<SiTrackerGSMatchedRecHit2DCollection> theGSRecHits;
00181 e.getByLabel(hitProducer, theGSRecHits);
00182
00183
00184 const std::vector<unsigned> theSimTrackIds = theGSRecHits->ids();
00185
00186
00187 edm::Handle<edm::SimVertexContainer> theSimVtx;
00188 e.getByLabel("famosSimHits",theSimVtx);
00189 edm::Handle<edm::SimTrackContainer> theSTC;
00190 e.getByLabel("famosSimHits",theSTC);
00191
00192 const edm::SimTrackContainer* theSimTracks = &(*theSTC);
00193 LogDebug("FastTracking")<<"looking at: "<< theSimTrackIds.size()<<" simtracks.";
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203
00204
00205 std::vector<edm::Handle<reco::TrackCollection> > theTrackCollections;
00206 std::vector<edm:: Handle<std::vector<Trajectory> > > theTrajectoryCollections;
00207 std::vector<edm::Handle<TrajTrackAssociationCollection> > theAssoMaps;
00208 std::vector<bool> isTrackCollections;
00209 TrajTrackAssociationCollection::const_iterator anAssociation;
00210 TrajTrackAssociationCollection::const_iterator lastAssociation;
00211 TrackMap theTrackMap;
00212 unsigned nCollections = trackProducers.size();
00213 unsigned nRecoHits = 0;
00214
00215 if ( nCollections ) {
00216 theTrackCollections.resize(nCollections);
00217 theTrajectoryCollections.resize(nCollections);
00218 theAssoMaps.resize(nCollections);
00219 isTrackCollections.resize(nCollections);
00220 for ( unsigned tprod=0; tprod < nCollections; ++tprod ) {
00221 isTrackCollections[tprod] = e.getByLabel(trackProducers[tprod],theTrackCollections[tprod]);
00222
00223 if ( isTrackCollections[tprod] ) {
00224
00225 reco::TrackCollection::const_iterator aTrack = theTrackCollections[tprod]->begin();
00226 reco::TrackCollection::const_iterator lastTrack = theTrackCollections[tprod]->end();
00227
00228 for ( ; aTrack!=lastTrack; ++aTrack ) nRecoHits+= aTrack->recHitsSize();
00229 e.getByLabel(trackProducers[tprod],theTrajectoryCollections[tprod]);
00230 e.getByLabel(trackProducers[tprod],theAssoMaps[tprod]);
00231
00232 anAssociation = theAssoMaps[tprod]->begin();
00233 lastAssociation = theAssoMaps[tprod]->end();
00234 #ifdef FAMOS_DEBUG
00235 std::cout << "Input Track Producer " << tprod << " : " << trackProducers[tprod] << std::endl;
00236 std::cout << "List of tracks already reconstructed " << std::endl;
00237 #endif
00238
00239 for ( ; anAssociation != lastAssociation; ++anAssociation ) {
00240 edm::Ref<std::vector<Trajectory> > aTrajectoryRef = anAssociation->key;
00241 reco::TrackRef aTrackRef = anAssociation->val;
00242
00243 int recoTrackId = findId(*aTrackRef);
00244 if ( recoTrackId < 0 ) continue;
00245 #ifdef FAMOS_DEBUG
00246 std::cout << recoTrackId << " ";
00247 #endif
00248
00249 theTrackMap[recoTrackId] = TrackPair(aTrackRef,aTrajectoryRef);
00250 }
00251 #ifdef FAMOS_DEBUG
00252 std::cout << std::endl;
00253 #endif
00254 }
00255 }
00256
00257 recoHits->reserve(nRecoHits);
00258 }
00259
00260
00261
00262 int currentTrackId = -1;
00263
00264
00265
00266
00267
00268
00269
00270 #ifdef FAMOS_DEBUG
00271 std::cout << "Input seed Producer : " << seedProducer << std::endl;
00272 std::cout << "Number of seeds : " << theSeeds->size() << std::endl;
00273 #endif
00274 unsigned seed_size = theSeeds->size();
00275 for (unsigned seednr = 0; seednr < seed_size; ++seednr){
00276
00277 LogDebug("FastTracking")<<"looking at seed #:"<<seednr;
00278
00279
00280 const BasicTrajectorySeed* aSeed = &((*theSeeds)[seednr]);
00281
00282 std::vector<int> simTrackIds;
00283 std::map<int,TrajectoryStateOnSurface> seedStates;
00284 std::map<int,TrajectoryStateOnSurface> simtkStates;
00285
00286 TrackerRecHit theFirstSeedingTrackerRecHit;
00287 if (theSeeds->at(seednr).nHits()==0){
00288
00289
00290 LogDebug("FastTracking")<<" seed with no hits to be considered.";
00291
00292
00293
00294
00295
00296 PTrajectoryStateOnDet ptod =theSeeds->at(seednr).startingState();
00297 DetId id(ptod.detId());
00298 const GeomDet * g = theGeometry->idToDet(id);
00299 const Surface * surface=&g->surface();
00300
00301 TrajectoryStateOnSurface seedState(trajectoryStateTransform::transientState(ptod,surface,theMagField));
00302
00303 edm::ESHandle<Propagator> propagator;
00304 es.get<TrackingComponentsRecord>().get("AnyDirectionAnalyticalPropagator",propagator);
00305
00306
00307
00308
00309
00310
00311
00312 double minimunEst=1000000;
00313 LogDebug("FastTracking")<<"looking at: "<< theSimTrackIds.size()<<" simtracks.";
00314 for ( unsigned tkId=0; tkId != theSimTrackIds.size(); ++tkId ) {
00315
00316 const SimTrack & simtrack = (*theSimTracks)[theSimTrackIds[tkId]];
00317
00318 GlobalPoint position(simtrack.trackerSurfacePosition().x(),
00319 simtrack.trackerSurfacePosition().y(),
00320 simtrack.trackerSurfacePosition().z());
00321
00322 GlobalVector momentum(simtrack.trackerSurfaceMomentum().x(),
00323 simtrack.trackerSurfaceMomentum().y(),
00324 simtrack.trackerSurfaceMomentum().z());
00325
00326 if (position.basicVector().dot( momentum.basicVector() ) * seedState.globalPosition().basicVector().dot(seedState.globalMomentum().basicVector()) <0. ){
00327 LogDebug("FastTracking")<<"not on the same direction.";
00328 continue;
00329 }
00330
00331
00332 int charge = (int) simtrack.charge();
00333 GlobalTrajectoryParameters glb_parameters(position,
00334 momentum,
00335 charge,
00336 theMagField);
00337 FreeTrajectoryState simtrack_trackerstate(glb_parameters);
00338
00339 TrajectoryStateOnSurface simtrack_comparestate = propagator->propagate(simtrack_trackerstate,*surface);
00340
00341
00342 if (!simtrack_comparestate.isValid()){
00343 LogDebug("FastTracking")<<" ok this is a state-based seed. simtrack state does not propagate to the seed surface. skipping.";
00344 continue;}
00345
00346 if (simtrack_comparestate.globalPosition().basicVector().dot(simtrack_comparestate.globalMomentum().basicVector()) * seedState.globalPosition().basicVector().dot(seedState.globalMomentum().basicVector()) <0. ){
00347 LogDebug("FastTracking")<<"not on the same direction.";
00348 continue;
00349 }
00350
00351
00352 AlgebraicVector5 v(seedState.localParameters().vector() - simtrack_comparestate.localParameters().vector());
00353 AlgebraicSymMatrix55 m(seedState.localError().matrix());
00354 bool ierr = !m.Invert();
00355 if ( ierr ){
00356 edm::LogWarning("FastTracking") <<" Candidate Producer cannot invert the error matrix! - Skipping...";
00357 continue;
00358 }
00359 double est = ROOT::Math::Similarity(v,m);
00360 LogDebug("FastTracking")<<"comparing two state on the seed surface:\n"
00361 <<"seed: "<<seedState
00362 <<"sim: "<<simtrack_comparestate
00363 <<"\n estimator is: "<<est;
00364
00365 if (est<minimunEst) minimunEst=est;
00366 if (est<estimatorCut_){
00367 simTrackIds.push_back(theSimTrackIds[tkId]);
00368
00369
00370
00371
00372 AlgebraicSymMatrix55 C = seedState.curvilinearError().matrix();
00373 C *= 0.0000001;
00374
00375 seedStates[theSimTrackIds[tkId]] = TrajectoryStateOnSurface(simtrack_comparestate.globalParameters(),
00376 CurvilinearTrajectoryError(C),
00377 seedState.surface());
00378 LogDebug("FastTracking")<<"the compatibility estimator is: "<<est<<" for track id: "<<simTrackIds.back();
00379 }
00380 }
00381 if (simTrackIds.size()==0) LogDebug("FastTracking")<<"could not find any simtrack within errors, closest was at: "<<minimunEst;
00382 }
00383 else{
00384
00385
00386 TrajectorySeed::range theSeedingRecHitRange = aSeed->recHits();
00387 const SiTrackerGSMatchedRecHit2D * theFirstSeedingRecHit = (const SiTrackerGSMatchedRecHit2D*) (&(*(theSeedingRecHitRange.first)));
00388 theFirstSeedingTrackerRecHit = TrackerRecHit(theFirstSeedingRecHit,theGeometry,tTopo);
00389
00390 simTrackIds.push_back( theFirstSeedingRecHit->simtrackId() );
00391 }
00392
00393
00394
00395
00396 for (unsigned int iToMake=0;iToMake!=simTrackIds.size();++iToMake){
00397 int simTrackId = simTrackIds[iToMake];
00398
00399
00400
00401 if ( seedCleaning && currentTrackId == simTrackId ) continue;
00402 currentTrackId = simTrackId;
00403
00404
00405 std::vector<TrackerRecHit> theTrackerRecHits;
00406 unsigned theNumberOfCrossedLayers = 0;
00407
00408
00409 TrackMap::const_iterator theTrackIt = theTrackMap.find(simTrackId);
00410 if ( nCollections && theTrackIt != theTrackMap.end() ) {
00411
00412 if ( keepFittedTracks ) {
00413 LogDebug("FastTracking") << "Track " << simTrackId << " already reconstructed -> copy it";
00414
00415 reco::TrackRef aTrackRef = theTrackIt->second.first;
00416 edm::Ref<std::vector<Trajectory> > aTrajectoryRef = theTrackIt->second.second;
00417
00418
00419 reco::Track aRecoTrack(*aTrackRef);
00420 recoTracks->push_back(aRecoTrack);
00421
00422
00423 unsigned nh = aRecoTrack.recHitsSize();
00424 for ( unsigned ih=0; ih<nh; ++ih ) {
00425 TrackingRecHit *hit = aRecoTrack.recHit(ih)->clone();
00426 recoHits->push_back(hit);
00427 }
00428
00429
00430 recoTrajectories->push_back(*aTrajectoryRef);
00431
00432 }
00433 else {
00434 LogDebug("FastTracking") << "Track " << simTrackId << " already reconstructed -> ignore it";
00435 }
00436
00437
00438
00439 }
00440 else{
00441
00442 LogDebug("FastTracking")<<"Track " << simTrackId << " is considered to return a track candidate" ;
00443
00444
00445 SiTrackerGSMatchedRecHit2DCollection::range theRecHitRange = theGSRecHits->get(simTrackId);
00446 SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorBegin = theRecHitRange.first;
00447 SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorEnd = theRecHitRange.second;
00448 SiTrackerGSMatchedRecHit2DCollection::const_iterator iterRecHit;
00449
00450 LogDebug("FastTracking")<<"counting: "<<theRecHitRangeIteratorEnd-theRecHitRangeIteratorBegin<<" hits to be considered.";
00451
00452 bool firstRecHit = true;
00453 TrackerRecHit theCurrentRecHit, thePreviousRecHit;
00454 TrackerRecHit theFirstHitComp, theSecondHitComp;
00455
00456 for ( iterRecHit = theRecHitRangeIteratorBegin;
00457 iterRecHit != theRecHitRangeIteratorEnd;
00458 ++iterRecHit) {
00459
00460
00461 if ( theNumberOfCrossedLayers >= maxNumberOfCrossedLayers ) continue;
00462
00463
00464 if(!firstRecHit) thePreviousRecHit = theCurrentRecHit;
00465 theCurrentRecHit = TrackerRecHit(&(*iterRecHit),theGeometry,tTopo);
00466
00467
00468
00469
00470
00471
00472
00473
00474
00475
00476
00477 if ( !theCurrentRecHit.isOnTheSameLayer(thePreviousRecHit) )
00478 ++theNumberOfCrossedLayers;
00479
00480
00481
00482 if ( !rejectOverlaps || firstRecHit ) {
00483
00484 if ( splitHits && theCurrentRecHit.matchedHit()->isMatched() ) addSplitHits(theCurrentRecHit,theTrackerRecHits);
00485 else theTrackerRecHits.push_back(theCurrentRecHit);
00486 firstRecHit = false;
00487
00488
00489
00490 } else {
00491
00492
00493 if ( theCurrentRecHit.subDetId() != thePreviousRecHit.subDetId() ||
00494 theCurrentRecHit.layerNumber() != thePreviousRecHit.layerNumber() ) {
00495
00496
00497 if ( splitHits && theCurrentRecHit.matchedHit()->isMatched() ) addSplitHits(theCurrentRecHit,theTrackerRecHits);
00498 else theTrackerRecHits.push_back(theCurrentRecHit);
00499
00500
00501 } else if ( theCurrentRecHit.localError() < thePreviousRecHit.localError() ) {
00502
00503
00504 if( splitHits && theCurrentRecHit.matchedHit()->isMatched() ){
00505
00506 theTrackerRecHits.pop_back();
00507 if ( thePreviousRecHit.matchedHit()->isMatched() ) theTrackerRecHits.pop_back();
00508
00509
00510 addSplitHits(theCurrentRecHit,theTrackerRecHits);
00511 }
00512
00513 else {
00514 theTrackerRecHits.back() = theCurrentRecHit;
00515 }
00516
00517 } else {
00518
00519
00520 theCurrentRecHit = thePreviousRecHit;
00521 }
00522 }
00523 }
00524
00525
00526 }
00527
00528 LogDebug("FastTracking")<<" number of hits: " << theTrackerRecHits.size()<<" after counting overlaps and splitting.";
00529
00530
00531 edm::OwnVector<TrackingRecHit> recHits;
00532 unsigned nTrackerHits = theTrackerRecHits.size();
00533 recHits.reserve(nTrackerHits);
00534
00535 if (aSeed->direction()==oppositeToMomentum){
00536 LogDebug("FastTracking")<<"reversing the order of the hits";
00537 std::reverse(theTrackerRecHits.begin(),theTrackerRecHits.end());
00538 }
00539
00540 for ( unsigned ih=0; ih<nTrackerHits; ++ih ) {
00541 TrackingRecHit* aTrackingRecHit = theTrackerRecHits[ih].hit()->clone();
00542 recHits.push_back(aTrackingRecHit);
00543
00544 const DetId& detId = theTrackerRecHits[ih].hit()->geographicalId();
00545 LogDebug("FastTracking")
00546 << "Added RecHit from detid " << detId.rawId()
00547 << " subdet = " << theTrackerRecHits[ih].subDetId()
00548 << " layer = " << theTrackerRecHits[ih].layerNumber()
00549 << " ring = " << theTrackerRecHits[ih].ringNumber()
00550 << " error = " << theTrackerRecHits[ih].localError()
00551 << std::endl
00552
00553 << "Track/z/r : "
00554 << simTrackId << " "
00555 << theTrackerRecHits[ih].globalPosition().z() << " "
00556 << theTrackerRecHits[ih].globalPosition().perp() << std::endl;
00557 if ( theTrackerRecHits[ih].matchedHit() && theTrackerRecHits[ih].matchedHit()->isMatched() )
00558 LogTrace("FastTracking") << "Matched : " << theTrackerRecHits[ih].matchedHit()->isMatched()
00559 << "Rphi Hit = " << theTrackerRecHits[ih].matchedHit()->monoHit()->simhitId()
00560 << "Stereo Hit = " << theTrackerRecHits[ih].matchedHit()->stereoHit()->simhitId()
00561 <<std::endl;
00562 }
00563
00564
00565 if ( theNumberOfCrossedLayers < minNumberOfCrossedLayers ) {
00566 LogDebug("FastTracking")<<"not enough layer crossed ("<<theNumberOfCrossedLayers<<")";
00567 continue;
00568 }
00569
00570
00571
00572
00573
00574
00575 PTrajectoryStateOnDet PTSOD;
00576
00577 if (aSeed->nHits()==0){
00578
00579
00580
00581 PTSOD = trajectoryStateTransform::persistentState(seedStates[simTrackId],aSeed->startingState().detId());
00582
00583 } else {
00584
00585 int vertexIndex = (*theSimTracks)[currentTrackId].vertIndex();
00586
00587 GlobalPoint position((*theSimVtx)[vertexIndex].position().x(),
00588 (*theSimVtx)[vertexIndex].position().y(),
00589 (*theSimVtx)[vertexIndex].position().z());
00590
00591
00592 GlobalVector momentum( (*theSimTracks)[currentTrackId].momentum().x() ,
00593 (*theSimTracks)[currentTrackId].momentum().y() ,
00594 (*theSimTracks)[currentTrackId].momentum().z() );
00595
00596 float charge = (*theSimTracks)[simTrackId].charge();
00597
00598 GlobalTrajectoryParameters initialParams(position,momentum,(int)charge,theMagField);
00599
00600 AlgebraicSymMatrix55 errorMatrix= AlgebraicMatrixID();
00601 CurvilinearTrajectoryError initialError(errorMatrix);
00602
00603 FreeTrajectoryState initialFTS(initialParams, initialError);
00604 #ifdef FAMOS_DEBUG
00605 std::cout << "TrajectorySeedProducer: FTS momentum " << initialFTS.momentum() << std::endl;
00606 #endif
00607 const GeomDet* initialLayer = theGeometry->idToDet(recHits.front().geographicalId());
00608
00609
00610 const TrajectoryStateOnSurface initialTSOS = thePropagator->propagate(initialFTS,initialLayer->surface()) ;
00611 if (!initialTSOS.isValid()) continue;
00612
00613
00614 PTSOD = trajectoryStateTransform::persistentState(initialTSOS,recHits.front().geographicalId().rawId());
00615 }
00616
00617 TrackCandidate
00618 newTrackCandidate(recHits,
00619 *aSeed,
00620 PTSOD,
00621 edm::RefToBase<TrajectorySeed>(theSeeds,seednr));
00622
00623 LogDebug("FastTracking")<< "\tSeed Information " << std::endl
00624 << "\tSeed Direction = " << aSeed->direction() << std::endl
00625 << "\tSeed StartingDet = " << aSeed->startingState().detId() << std::endl
00626 << "\tTrajectory Parameters " << std::endl
00627 << "\t\t detId = " << newTrackCandidate.trajectoryStateOnDet().detId() << std::endl
00628 << "\t\t loc.px = "
00629 << newTrackCandidate.trajectoryStateOnDet().parameters().momentum().x()
00630 << std::endl
00631 << "\t\t loc.py = "
00632 << newTrackCandidate.trajectoryStateOnDet().parameters().momentum().y()
00633 << std::endl
00634 << "\t\t loc.pz = "
00635 << newTrackCandidate.trajectoryStateOnDet().parameters().momentum().z()
00636 << std::endl
00637 << "\t\t error = ";
00638
00639 bool newTrackCandidateIsDuplicate = isDuplicateCandidate(*output,newTrackCandidate);
00640 if (!newTrackCandidateIsDuplicate) output->push_back(newTrackCandidate);
00641 LogDebug("FastTracking")<<"filling a track candidate into the collection, now having: "<<output->size();
00642
00643 }
00644 }
00645
00646
00647 LogDebug("FastTracking") << "Saving "
00648 << output->size() << " track candidates and "
00649 << recoTracks->size() << " reco::Tracks ";
00650
00651 e.put(output);
00652
00653
00654
00655
00656
00657 edm::OrphanHandle <TrackingRecHitCollection> theRecoHits = e.put(recoHits );
00658
00659
00660 unsigned hits=0;
00661 unsigned nTracks = recoTracks->size();
00662 recoTrackExtras->reserve(nTracks);
00663 for ( unsigned index = 0; index < nTracks; ++index ) {
00664
00665 reco::Track& aTrack = recoTracks->at(index);
00666 reco::TrackExtra aTrackExtra(aTrack.outerPosition(),
00667 aTrack.outerMomentum(),
00668 aTrack.outerOk(),
00669 aTrack.innerPosition(),
00670 aTrack.innerMomentum(),
00671 aTrack.innerOk(),
00672 aTrack.outerStateCovariance(),
00673 aTrack.outerDetId(),
00674 aTrack.innerStateCovariance(),
00675 aTrack.innerDetId(),
00676 aTrack.seedDirection(),
00677 aTrack.seedRef());
00678
00679 unsigned nHits = aTrack.recHitsSize();
00680 for ( unsigned int ih=0; ih<nHits; ++ih) {
00681 aTrackExtra.add(TrackingRecHitRef(theRecoHits,hits++));
00682 }
00683 recoTrackExtras->push_back(aTrackExtra);
00684 }
00685
00686
00687
00688 edm::OrphanHandle<reco::TrackExtraCollection> theRecoTrackExtras = e.put(recoTrackExtras);
00689
00690
00691 for ( unsigned index = 0; index<nTracks; ++index ) {
00692 const reco::TrackExtraRef theTrackExtraRef(theRecoTrackExtras,index);
00693 (recoTracks->at(index)).setExtra(theTrackExtraRef);
00694 }
00695
00696
00697 edm::OrphanHandle<reco::TrackCollection> theRecoTracks = e.put(recoTracks);
00698
00699
00700 edm::OrphanHandle<std::vector<Trajectory> > theRecoTrajectories = e.put(recoTrajectories);
00701
00702
00703 for ( unsigned index = 0; index<nTracks; ++index ) {
00704 edm::Ref<std::vector<Trajectory> > trajRef( theRecoTrajectories, index );
00705 edm::Ref<reco::TrackCollection> tkRef( theRecoTracks, index );
00706 recoTrajTrackMap->insert(trajRef,tkRef);
00707 }
00708
00709
00710
00711 e.put(recoTrajTrackMap);
00712
00713 }
00714
00715 int
00716 TrackCandidateProducer::findId(const reco::Track& aTrack) const {
00717 int trackId = -1;
00718 trackingRecHit_iterator aHit = aTrack.recHitsBegin();
00719 trackingRecHit_iterator lastHit = aTrack.recHitsEnd();
00720 for ( ; aHit!=lastHit; ++aHit ) {
00721 if ( !aHit->get()->isValid() ) continue;
00722
00723 const SiTrackerGSMatchedRecHit2D * rechit = (const SiTrackerGSMatchedRecHit2D*) (aHit->get());
00724 trackId = rechit->simtrackId();
00725 break;
00726 }
00727 return trackId;
00728 }
00729
00730 void
00731 TrackCandidateProducer::addSplitHits(const TrackerRecHit& theCurrentRecHit,
00732 std::vector<TrackerRecHit>& theTrackerRecHits) {
00733
00734 const SiTrackerGSRecHit2D* mHit = theCurrentRecHit.matchedHit()->monoHit();
00735 const SiTrackerGSRecHit2D* sHit = theCurrentRecHit.matchedHit()->stereoHit();
00736
00737
00738 if( mHit->simhitId() < sHit->simhitId() ) {
00739
00740 theTrackerRecHits.push_back(TrackerRecHit(mHit,theCurrentRecHit));
00741 theTrackerRecHits.push_back(TrackerRecHit(sHit,theCurrentRecHit));
00742
00743 } else {
00744
00745 theTrackerRecHits.push_back(TrackerRecHit(sHit,theCurrentRecHit));
00746 theTrackerRecHits.push_back(TrackerRecHit(mHit,theCurrentRecHit));
00747
00748 }
00749
00750 }
00751
00752 bool TrackCandidateProducer::isDuplicateCandidate(const TrackCandidateCollection& candidates, const TrackCandidate& newCand) const{
00753 typedef TrackCandidateCollection::const_iterator TCCI;
00754
00755 TCCI candsEnd = candidates.end();
00756 double newQbp = newCand.trajectoryStateOnDet().parameters().qbp();
00757 double newDxdz = newCand.trajectoryStateOnDet().parameters().dxdz();
00758 TrackCandidate::range newHits = newCand.recHits();
00759
00760 for (TCCI iCand = candidates.begin(); iCand!= candsEnd; ++iCand){
00761
00762 double iQbp = iCand->trajectoryStateOnDet().parameters().qbp();
00763 double iDxdz = iCand->trajectoryStateOnDet().parameters().dxdz();
00764 if (newQbp == iQbp && newDxdz == iDxdz){
00765 LogDebug("isDuplicateCandidate")<<"Probably a duplicate "<<iQbp <<" : "<<iDxdz;
00766 TrackCandidate::range iHits = iCand->recHits();
00767 unsigned int nHits = 0;
00768 unsigned int nShared = 0;
00769 unsigned int nnewHits = 0;
00770 for (TrackCandidate::const_iterator iiHit = iHits.first; iiHit != iHits.second; ++iiHit){
00771 nHits++;
00772 for (TrackCandidate::const_iterator inewHit = newHits.first; inewHit != newHits.second; ++inewHit){
00773 if (iiHit == iHits.first) nnewHits++;
00774 if (sameLocalParameters(&*iiHit,&*inewHit)) nShared++;
00775 }
00776 }
00777 LogDebug("isDuplicateCandidate") <<"nHits "<<nHits<<" nnewHits "<< nnewHits<< " nShared "<<nShared;
00778 if (nHits == nShared && nHits == nnewHits) return true;
00779 }
00780 }
00781 return false;
00782 }
00783
00784 bool TrackCandidateProducer::sameLocalParameters(const TrackingRecHit* aH, const TrackingRecHit* bH) const {
00785 bool localPointSame = aH->localPosition() == bH->localPosition();
00786 bool localErrorSame = aH->localPositionError().xx() == bH->localPositionError().xx();
00787 localErrorSame &= aH->localPositionError().yy() == bH->localPositionError().yy();
00788 localErrorSame &= aH->localPositionError().xy() == bH->localPositionError().xy();
00789 return localPointSame && localErrorSame;
00790 }