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 & run, 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 if(theSeeds->size() == 0) {
00163 e.put(output);
00164 e.put(recoTracks);
00165 e.put(recoHits);
00166 e.put(recoTrackExtras);
00167 e.put(recoTrajectories);
00168 e.put(recoTrajTrackMap);
00169 return;
00170 }
00171
00172
00173
00174 edm::Handle<SiTrackerGSMatchedRecHit2DCollection> theGSRecHits;
00175 e.getByLabel(hitProducer, theGSRecHits);
00176
00177
00178 const std::vector<unsigned> theSimTrackIds = theGSRecHits->ids();
00179
00180
00181 edm::Handle<edm::SimVertexContainer> theSimVtx;
00182 e.getByLabel("famosSimHits",theSimVtx);
00183 edm::Handle<edm::SimTrackContainer> theSTC;
00184 e.getByLabel("famosSimHits",theSTC);
00185
00186 const edm::SimTrackContainer* theSimTracks = &(*theSTC);
00187 LogDebug("FastTracking")<<"looking at: "<< theSimTrackIds.size()<<" simtracks.";
00188
00189
00190
00191
00192
00193
00194
00195
00196
00197
00198
00199 std::vector<edm::Handle<reco::TrackCollection> > theTrackCollections;
00200 std::vector<edm:: Handle<std::vector<Trajectory> > > theTrajectoryCollections;
00201 std::vector<edm::Handle<TrajTrackAssociationCollection> > theAssoMaps;
00202 std::vector<bool> isTrackCollections;
00203 TrajTrackAssociationCollection::const_iterator anAssociation;
00204 TrajTrackAssociationCollection::const_iterator lastAssociation;
00205 TrackMap theTrackMap;
00206 unsigned nCollections = trackProducers.size();
00207 unsigned nRecoHits = 0;
00208
00209 if ( nCollections ) {
00210 theTrackCollections.resize(nCollections);
00211 theTrajectoryCollections.resize(nCollections);
00212 theAssoMaps.resize(nCollections);
00213 isTrackCollections.resize(nCollections);
00214 for ( unsigned tprod=0; tprod < nCollections; ++tprod ) {
00215 isTrackCollections[tprod] = e.getByLabel(trackProducers[tprod],theTrackCollections[tprod]);
00216
00217 if ( isTrackCollections[tprod] ) {
00218
00219 reco::TrackCollection::const_iterator aTrack = theTrackCollections[tprod]->begin();
00220 reco::TrackCollection::const_iterator lastTrack = theTrackCollections[tprod]->end();
00221
00222 for ( ; aTrack!=lastTrack; ++aTrack ) nRecoHits+= aTrack->recHitsSize();
00223 e.getByLabel(trackProducers[tprod],theTrajectoryCollections[tprod]);
00224 e.getByLabel(trackProducers[tprod],theAssoMaps[tprod]);
00225
00226 anAssociation = theAssoMaps[tprod]->begin();
00227 lastAssociation = theAssoMaps[tprod]->end();
00228 #ifdef FAMOS_DEBUG
00229 std::cout << "Input Track Producer " << tprod << " : " << trackProducers[tprod] << std::endl;
00230 std::cout << "List of tracks already reconstructed " << std::endl;
00231 #endif
00232
00233 for ( ; anAssociation != lastAssociation; ++anAssociation ) {
00234 edm::Ref<std::vector<Trajectory> > aTrajectoryRef = anAssociation->key;
00235 reco::TrackRef aTrackRef = anAssociation->val;
00236
00237 int recoTrackId = findId(*aTrackRef);
00238 if ( recoTrackId < 0 ) continue;
00239 #ifdef FAMOS_DEBUG
00240 std::cout << recoTrackId << " ";
00241 #endif
00242
00243 theTrackMap[recoTrackId] = TrackPair(aTrackRef,aTrajectoryRef);
00244 }
00245 #ifdef FAMOS_DEBUG
00246 std::cout << std::endl;
00247 #endif
00248 }
00249 }
00250
00251 recoHits->reserve(nRecoHits);
00252 }
00253
00254
00255
00256 int currentTrackId = -1;
00257
00258
00259
00260
00261
00262
00263
00264 #ifdef FAMOS_DEBUG
00265 std::cout << "Input seed Producer : " << seedProducer << std::endl;
00266 std::cout << "Number of seeds : " << theSeeds->size() << std::endl;
00267 #endif
00268 unsigned seed_size = theSeeds->size();
00269 for (unsigned seednr = 0; seednr < seed_size; ++seednr){
00270
00271 LogDebug("FastTracking")<<"looking at seed #:"<<seednr;
00272
00273
00274 const BasicTrajectorySeed* aSeed = &((*theSeeds)[seednr]);
00275
00276 std::vector<int> simTrackIds;
00277 std::map<int,TrajectoryStateOnSurface> seedStates;
00278 std::map<int,TrajectoryStateOnSurface> simtkStates;
00279
00280 TrackerRecHit theFirstSeedingTrackerRecHit;
00281 if (theSeeds->at(seednr).nHits()==0){
00282
00283
00284 LogDebug("FastTracking")<<" seed with no hits to be considered.";
00285
00286
00287
00288
00289
00290 PTrajectoryStateOnDet ptod =theSeeds->at(seednr).startingState();
00291 DetId id(ptod.detId());
00292 const GeomDet * g = theGeometry->idToDet(id);
00293 const Surface * surface=&g->surface();
00294 TrajectoryStateTransform tsTransform;
00295 TrajectoryStateOnSurface seedState(tsTransform.transientState(ptod,surface,theMagField));
00296
00297 edm::ESHandle<Propagator> propagator;
00298 es.get<TrackingComponentsRecord>().get("AnyDirectionAnalyticalPropagator",propagator);
00299
00300
00301
00302
00303
00304
00305
00306 double minimunEst=1000000;
00307 LogDebug("FastTracking")<<"looking at: "<< theSimTrackIds.size()<<" simtracks.";
00308 for ( unsigned tkId=0; tkId != theSimTrackIds.size(); ++tkId ) {
00309
00310 const SimTrack & simtrack = (*theSimTracks)[theSimTrackIds[tkId]];
00311
00312 GlobalPoint position(simtrack.trackerSurfacePosition().x(),
00313 simtrack.trackerSurfacePosition().y(),
00314 simtrack.trackerSurfacePosition().z());
00315
00316 GlobalVector momentum(simtrack.trackerSurfaceMomentum().x(),
00317 simtrack.trackerSurfaceMomentum().y(),
00318 simtrack.trackerSurfaceMomentum().z());
00319
00320 if (position.basicVector().dot( momentum.basicVector() ) * seedState.globalPosition().basicVector().dot(seedState.globalMomentum().basicVector()) <0. ){
00321 LogDebug("FastTracking")<<"not on the same direction.";
00322 continue;
00323 }
00324
00325
00326 int charge = (int) simtrack.charge();
00327 GlobalTrajectoryParameters glb_parameters(position,
00328 momentum,
00329 charge,
00330 theMagField);
00331 FreeTrajectoryState simtrack_trackerstate(glb_parameters);
00332
00333 TrajectoryStateOnSurface simtrack_comparestate = propagator->propagate(simtrack_trackerstate,*surface);
00334
00335
00336 if (!simtrack_comparestate.isValid()){
00337 LogDebug("FastTracking")<<" ok this is a state-based seed. simtrack state does not propagate to the seed surface. skipping.";
00338 continue;}
00339
00340 if (simtrack_comparestate.globalPosition().basicVector().dot(simtrack_comparestate.globalMomentum().basicVector()) * seedState.globalPosition().basicVector().dot(seedState.globalMomentum().basicVector()) <0. ){
00341 LogDebug("FastTracking")<<"not on the same direction.";
00342 continue;
00343 }
00344
00345
00346 AlgebraicVector5 v(seedState.localParameters().vector() - simtrack_comparestate.localParameters().vector());
00347 AlgebraicSymMatrix55 m(seedState.localError().matrix());
00348 bool ierr = !m.Invert();
00349 if ( ierr ){
00350 edm::LogWarning("FastTracking") <<" Candidate Producer cannot invert the error matrix! - Skipping...";
00351 continue;
00352 }
00353 double est = ROOT::Math::Similarity(v,m);
00354 LogDebug("FastTracking")<<"comparing two state on the seed surface:\n"
00355 <<"seed: "<<seedState
00356 <<"sim: "<<simtrack_comparestate
00357 <<"\n estimator is: "<<est;
00358
00359 if (est<minimunEst) minimunEst=est;
00360 if (est<estimatorCut_){
00361 simTrackIds.push_back(theSimTrackIds[tkId]);
00362
00363
00364
00365
00366 AlgebraicSymMatrix55 C = seedState.curvilinearError().matrix();
00367 C *= 0.0000001;
00368
00369 seedStates[theSimTrackIds[tkId]] = TrajectoryStateOnSurface(simtrack_comparestate.globalParameters(),
00370 CurvilinearTrajectoryError(C),
00371 seedState.surface());
00372 LogDebug("FastTracking")<<"the compatibility estimator is: "<<est<<" for track id: "<<simTrackIds.back();
00373 }
00374 }
00375 if (simTrackIds.size()==0) LogDebug("FastTracking")<<"could not find any simtrack within errors, closest was at: "<<minimunEst;
00376 }
00377 else{
00378
00379
00380 TrajectorySeed::range theSeedingRecHitRange = aSeed->recHits();
00381 const SiTrackerGSMatchedRecHit2D * theFirstSeedingRecHit = (const SiTrackerGSMatchedRecHit2D*) (&(*(theSeedingRecHitRange.first)));
00382 theFirstSeedingTrackerRecHit = TrackerRecHit(theFirstSeedingRecHit,theGeometry);
00383
00384 simTrackIds.push_back( theFirstSeedingRecHit->simtrackId() );
00385 }
00386
00387
00388
00389
00390 for (unsigned int iToMake=0;iToMake!=simTrackIds.size();++iToMake){
00391 int simTrackId = simTrackIds[iToMake];
00392
00393
00394
00395 if ( seedCleaning && currentTrackId == simTrackId ) continue;
00396 currentTrackId = simTrackId;
00397
00398
00399 std::vector<TrackerRecHit> theTrackerRecHits;
00400 unsigned theNumberOfCrossedLayers = 0;
00401
00402
00403 TrackMap::const_iterator theTrackIt = theTrackMap.find(simTrackId);
00404 if ( nCollections && theTrackIt != theTrackMap.end() ) {
00405
00406 if ( keepFittedTracks ) {
00407 LogDebug("FastTracking") << "Track " << simTrackId << " already reconstructed -> copy it";
00408
00409 reco::TrackRef aTrackRef = theTrackIt->second.first;
00410 edm::Ref<std::vector<Trajectory> > aTrajectoryRef = theTrackIt->second.second;
00411
00412
00413 reco::Track aRecoTrack(*aTrackRef);
00414 recoTracks->push_back(aRecoTrack);
00415
00416
00417 unsigned nh = aRecoTrack.recHitsSize();
00418 for ( unsigned ih=0; ih<nh; ++ih ) {
00419 TrackingRecHit *hit = aRecoTrack.recHit(ih)->clone();
00420 recoHits->push_back(hit);
00421 }
00422
00423
00424 recoTrajectories->push_back(*aTrajectoryRef);
00425
00426 }
00427 else {
00428 LogDebug("FastTracking") << "Track " << simTrackId << " already reconstructed -> ignore it";
00429 }
00430
00431
00432
00433 }
00434 else{
00435
00436 LogDebug("FastTracking")<<"Track " << simTrackId << " is considered to return a track candidate" ;
00437
00438
00439 SiTrackerGSMatchedRecHit2DCollection::range theRecHitRange = theGSRecHits->get(simTrackId);
00440 SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorBegin = theRecHitRange.first;
00441 SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorEnd = theRecHitRange.second;
00442 SiTrackerGSMatchedRecHit2DCollection::const_iterator iterRecHit;
00443
00444 LogDebug("FastTracking")<<"counting: "<<theRecHitRangeIteratorEnd-theRecHitRangeIteratorBegin<<" hits to be considered.";
00445
00446 bool firstRecHit = true;
00447 TrackerRecHit theCurrentRecHit, thePreviousRecHit;
00448 TrackerRecHit theFirstHitComp, theSecondHitComp;
00449
00450 for ( iterRecHit = theRecHitRangeIteratorBegin;
00451 iterRecHit != theRecHitRangeIteratorEnd;
00452 ++iterRecHit) {
00453
00454
00455 if ( theNumberOfCrossedLayers >= maxNumberOfCrossedLayers ) continue;
00456
00457
00458 if(!firstRecHit) thePreviousRecHit = theCurrentRecHit;
00459 theCurrentRecHit = TrackerRecHit(&(*iterRecHit),theGeometry);
00460
00461
00462
00463
00464
00465
00466
00467
00468
00469
00470
00471 if ( !theCurrentRecHit.isOnTheSameLayer(thePreviousRecHit) )
00472 ++theNumberOfCrossedLayers;
00473
00474
00475
00476 if ( !rejectOverlaps || firstRecHit ) {
00477
00478 if ( splitHits && theCurrentRecHit.matchedHit()->isMatched() ) addSplitHits(theCurrentRecHit,theTrackerRecHits);
00479 else theTrackerRecHits.push_back(theCurrentRecHit);
00480 firstRecHit = false;
00481
00482
00483
00484 } else {
00485
00486
00487 if ( theCurrentRecHit.subDetId() != thePreviousRecHit.subDetId() ||
00488 theCurrentRecHit.layerNumber() != thePreviousRecHit.layerNumber() ) {
00489
00490
00491 if ( splitHits && theCurrentRecHit.matchedHit()->isMatched() ) addSplitHits(theCurrentRecHit,theTrackerRecHits);
00492 else theTrackerRecHits.push_back(theCurrentRecHit);
00493
00494
00495 } else if ( theCurrentRecHit.localError() < thePreviousRecHit.localError() ) {
00496
00497
00498 if( splitHits && theCurrentRecHit.matchedHit()->isMatched() ){
00499
00500 theTrackerRecHits.pop_back();
00501 if ( thePreviousRecHit.matchedHit()->isMatched() ) theTrackerRecHits.pop_back();
00502
00503
00504 addSplitHits(theCurrentRecHit,theTrackerRecHits);
00505 }
00506
00507 else {
00508 theTrackerRecHits.back() = theCurrentRecHit;
00509 }
00510
00511 } else {
00512
00513
00514 theCurrentRecHit = thePreviousRecHit;
00515 }
00516 }
00517 }
00518
00519
00520 }
00521
00522 LogDebug("FastTracking")<<" number of hits: " << theTrackerRecHits.size()<<" after counting overlaps and splitting.";
00523
00524
00525 edm::OwnVector<TrackingRecHit> recHits;
00526 unsigned nTrackerHits = theTrackerRecHits.size();
00527 recHits.reserve(nTrackerHits);
00528
00529 if (aSeed->direction()==oppositeToMomentum){
00530 LogDebug("FastTracking")<<"reversing the order of the hits";
00531 std::reverse(theTrackerRecHits.begin(),theTrackerRecHits.end());
00532 }
00533
00534 for ( unsigned ih=0; ih<nTrackerHits; ++ih ) {
00535 TrackingRecHit* aTrackingRecHit = theTrackerRecHits[ih].hit()->clone();
00536 recHits.push_back(aTrackingRecHit);
00537
00538 const DetId& detId = theTrackerRecHits[ih].hit()->geographicalId();
00539 LogDebug("FastTracking")
00540 << "Added RecHit from detid " << detId.rawId()
00541 << " subdet = " << theTrackerRecHits[ih].subDetId()
00542 << " layer = " << theTrackerRecHits[ih].layerNumber()
00543 << " ring = " << theTrackerRecHits[ih].ringNumber()
00544 << " error = " << theTrackerRecHits[ih].localError()
00545 << std::endl
00546
00547 << "Track/z/r : "
00548 << simTrackId << " "
00549 << theTrackerRecHits[ih].globalPosition().z() << " "
00550 << theTrackerRecHits[ih].globalPosition().perp() << std::endl;
00551 if ( theTrackerRecHits[ih].matchedHit() && theTrackerRecHits[ih].matchedHit()->isMatched() )
00552 LogTrace("FastTracking") << "Matched : " << theTrackerRecHits[ih].matchedHit()->isMatched()
00553 << "Rphi Hit = " << theTrackerRecHits[ih].matchedHit()->monoHit()->simhitId()
00554 << "Stereo Hit = " << theTrackerRecHits[ih].matchedHit()->stereoHit()->simhitId()
00555 <<std::endl;
00556 }
00557
00558
00559 if ( theNumberOfCrossedLayers < minNumberOfCrossedLayers ) {
00560 LogDebug("FastTracking")<<"not enough layer crossed ("<<theNumberOfCrossedLayers<<")";
00561 continue;
00562 }
00563
00564
00565
00566
00567
00568
00569 PTrajectoryStateOnDet PTSOD;
00570
00571 if (aSeed->nHits()==0){
00572
00573
00574 TrajectoryStateTransform tsTransform;
00575
00576 PTrajectoryStateOnDet * aPointer = tsTransform.persistentState(seedStates[simTrackId],aSeed->startingState().detId());
00577 PTSOD = *aPointer;
00578 delete aPointer;
00579
00580 } else {
00581
00582 int vertexIndex = (*theSimTracks)[currentTrackId].vertIndex();
00583
00584 GlobalPoint position((*theSimVtx)[vertexIndex].position().x(),
00585 (*theSimVtx)[vertexIndex].position().y(),
00586 (*theSimVtx)[vertexIndex].position().z());
00587
00588
00589 GlobalVector momentum( (*theSimTracks)[currentTrackId].momentum().x() ,
00590 (*theSimTracks)[currentTrackId].momentum().y() ,
00591 (*theSimTracks)[currentTrackId].momentum().z() );
00592
00593 float charge = (*theSimTracks)[simTrackId].charge();
00594
00595 GlobalTrajectoryParameters initialParams(position,momentum,(int)charge,theMagField);
00596
00597 AlgebraicSymMatrix55 errorMatrix= AlgebraicMatrixID();
00598 CurvilinearTrajectoryError initialError(errorMatrix);
00599
00600 FreeTrajectoryState initialFTS(initialParams, initialError);
00601 #ifdef FAMOS_DEBUG
00602 std::cout << "TrajectorySeedProducer: FTS momentum " << initialFTS.momentum() << std::endl;
00603 #endif
00604 const GeomDet* initialLayer = theGeometry->idToDet(recHits.front().geographicalId());
00605
00606
00607 const TrajectoryStateOnSurface initialTSOS = thePropagator->propagate(initialFTS,initialLayer->surface()) ;
00608 if (!initialTSOS.isValid()) continue;
00609 TrajectoryStateTransform tsTransform;
00610
00611 PTrajectoryStateOnDet * aPointer = tsTransform.persistentState(initialTSOS,recHits.front().geographicalId().rawId());
00612 PTSOD = *aPointer;
00613 delete aPointer;
00614 }
00615
00616 TrackCandidate
00617 newTrackCandidate(recHits,
00618 *aSeed,
00619 PTSOD,
00620 edm::RefToBase<TrajectorySeed>(theSeeds,seednr));
00621
00622 LogDebug("FastTracking")<< "\tSeed Information " << std::endl
00623 << "\tSeed Direction = " << aSeed->direction() << std::endl
00624 << "\tSeed StartingDet = " << aSeed->startingState().detId() << std::endl
00625 << "\tTrajectory Parameters " << std::endl
00626 << "\t\t detId = " << newTrackCandidate.trajectoryStateOnDet().detId() << std::endl
00627 << "\t\t loc.px = "
00628 << newTrackCandidate.trajectoryStateOnDet().parameters().momentum().x()
00629 << std::endl
00630 << "\t\t loc.py = "
00631 << newTrackCandidate.trajectoryStateOnDet().parameters().momentum().y()
00632 << std::endl
00633 << "\t\t loc.pz = "
00634 << newTrackCandidate.trajectoryStateOnDet().parameters().momentum().z()
00635 << std::endl
00636 << "\t\t error = ";
00637
00638
00639
00640
00641
00642
00643
00644 output->push_back(newTrackCandidate);
00645 LogDebug("FastTracking")<<"filling a track candidate into the collection, now having: "<<output->size();
00646
00647 }
00648 }
00649
00650
00651 LogDebug("FastTracking") << "Saving "
00652 << output->size() << " track candidates and "
00653 << recoTracks->size() << " reco::Tracks ";
00654
00655 e.put(output);
00656
00657
00658
00659
00660
00661 edm::OrphanHandle <TrackingRecHitCollection> theRecoHits = e.put(recoHits );
00662
00663
00664 unsigned hits=0;
00665 unsigned nTracks = recoTracks->size();
00666 recoTrackExtras->reserve(nTracks);
00667 for ( unsigned index = 0; index < nTracks; ++index ) {
00668
00669 reco::Track& aTrack = recoTracks->at(index);
00670 reco::TrackExtra aTrackExtra(aTrack.outerPosition(),
00671 aTrack.outerMomentum(),
00672 aTrack.outerOk(),
00673 aTrack.innerPosition(),
00674 aTrack.innerMomentum(),
00675 aTrack.innerOk(),
00676 aTrack.outerStateCovariance(),
00677 aTrack.outerDetId(),
00678 aTrack.innerStateCovariance(),
00679 aTrack.innerDetId(),
00680 aTrack.seedDirection(),
00681 aTrack.seedRef());
00682
00683 unsigned nHits = aTrack.recHitsSize();
00684 for ( unsigned int ih=0; ih<nHits; ++ih) {
00685 aTrackExtra.add(TrackingRecHitRef(theRecoHits,hits++));
00686 }
00687 recoTrackExtras->push_back(aTrackExtra);
00688 }
00689
00690
00691
00692 edm::OrphanHandle<reco::TrackExtraCollection> theRecoTrackExtras = e.put(recoTrackExtras);
00693
00694
00695 for ( unsigned index = 0; index<nTracks; ++index ) {
00696 const reco::TrackExtraRef theTrackExtraRef(theRecoTrackExtras,index);
00697 (recoTracks->at(index)).setExtra(theTrackExtraRef);
00698 }
00699
00700
00701 edm::OrphanHandle<reco::TrackCollection> theRecoTracks = e.put(recoTracks);
00702
00703
00704 edm::OrphanHandle<std::vector<Trajectory> > theRecoTrajectories = e.put(recoTrajectories);
00705
00706
00707 for ( unsigned index = 0; index<nTracks; ++index ) {
00708 edm::Ref<std::vector<Trajectory> > trajRef( theRecoTrajectories, index );
00709 edm::Ref<reco::TrackCollection> tkRef( theRecoTracks, index );
00710 recoTrajTrackMap->insert(trajRef,tkRef);
00711 }
00712
00713
00714
00715 e.put(recoTrajTrackMap);
00716
00717 }
00718
00719 int
00720 TrackCandidateProducer::findId(const reco::Track& aTrack) const {
00721 int trackId = -1;
00722 trackingRecHit_iterator aHit = aTrack.recHitsBegin();
00723 trackingRecHit_iterator lastHit = aTrack.recHitsEnd();
00724 for ( ; aHit!=lastHit; ++aHit ) {
00725 if ( !aHit->get()->isValid() ) continue;
00726
00727 const SiTrackerGSMatchedRecHit2D * rechit = (const SiTrackerGSMatchedRecHit2D*) (aHit->get());
00728 trackId = rechit->simtrackId();
00729 break;
00730 }
00731 return trackId;
00732 }
00733
00734 void
00735 TrackCandidateProducer::addSplitHits(const TrackerRecHit& theCurrentRecHit,
00736 std::vector<TrackerRecHit>& theTrackerRecHits) {
00737
00738 const SiTrackerGSRecHit2D* mHit = theCurrentRecHit.matchedHit()->monoHit();
00739 const SiTrackerGSRecHit2D* sHit = theCurrentRecHit.matchedHit()->stereoHit();
00740
00741
00742 if( mHit->simhitId() < sHit->simhitId() ) {
00743
00744 theTrackerRecHits.push_back(TrackerRecHit(mHit,theCurrentRecHit));
00745 theTrackerRecHits.push_back(TrackerRecHit(sHit,theCurrentRecHit));
00746
00747 } else {
00748
00749 theTrackerRecHits.push_back(TrackerRecHit(sHit,theCurrentRecHit));
00750 theTrackerRecHits.push_back(TrackerRecHit(mHit,theCurrentRecHit));
00751
00752 }
00753
00754 }