CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/RecoTracker/FinalTrackSelectors/src/DuplicateListMerger.cc

Go to the documentation of this file.
00001 #include "RecoTracker/FinalTrackSelectors/interface/DuplicateListMerger.h"
00002 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00003 #include "TrackingTools/PatternTools/interface/TrajTrackAssociation.h"
00004 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
00005 #include "RecoTracker/TrackProducer/interface/ClusterRemovalRefSetter.h"
00006 
00007 #include "FWCore/Framework/interface/Event.h"
00008 #include "CommonTools/Statistics/interface/ChiSquaredProbability.h"
00009 
00010 using reco::modules::DuplicateListMerger;
00011 
00012 DuplicateListMerger::DuplicateListMerger(const edm::ParameterSet& iPara) 
00013 {
00014   diffHitsCut_ = 9999;
00015   minTrkProbCut_ = 0.0;
00016   if(iPara.exists("diffHitsCut"))diffHitsCut_ = iPara.getParameter<int>("diffHitsCut");
00017   if(iPara.exists("minTrkProbCut"))minTrkProbCut_ = iPara.getParameter<double>("minTrkProbCut");
00018   if(iPara.exists("mergedSource"))mergedTrackSource_ = iPara.getParameter<edm::InputTag>("mergedSource");
00019   if(iPara.exists("originalSource"))originalTrackSource_ = iPara.getParameter<edm::InputTag>("originalSource");
00020   if(iPara.exists("candidateSource"))candidateSource_ = iPara.getParameter<edm::InputTag>("candidateSource");
00021   copyExtras_ = iPara.getUntrackedParameter<bool>("copyExtras",true);
00022   qualityToSet_ = reco::TrackBase::undefQuality;
00023   if (iPara.exists("newQuality")) {
00024     std::string qualityStr = iPara.getParameter<std::string>("newQuality");
00025     if (qualityStr != "") {
00026       qualityToSet_ = reco::TrackBase::qualityByName(qualityStr);
00027     }
00028   }
00029 
00030   produces<std::vector<reco::Track> >();
00031   produces< std::vector<Trajectory> >();
00032   produces< TrajTrackAssociationCollection >();
00033   makeReKeyedSeeds_ = iPara.getUntrackedParameter<bool>("makeReKeyedSeeds",false);
00034   if (makeReKeyedSeeds_){
00035     copyExtras_=true;
00036     produces<TrajectorySeedCollection>();
00037   }
00038   if(copyExtras_){
00039     produces<reco::TrackExtraCollection>();
00040     produces<TrackingRecHitCollection>();
00041   }
00042 
00043 }
00044 
00045 DuplicateListMerger::~DuplicateListMerger()
00046 {
00047 
00048   /* no op */
00049 
00050 }
00051 
00052 void DuplicateListMerger::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00053 {
00054   edm::Handle<edm::View<reco::Track> > originalHandle;
00055   iEvent.getByLabel(originalTrackSource_,originalHandle);
00056   edm::Handle<edm::View<reco::Track> > mergedHandle;
00057   iEvent.getByLabel(mergedTrackSource_,mergedHandle);
00058 
00059   edm::Handle< std::vector<Trajectory> >  mergedTrajHandle;
00060   iEvent.getByLabel(mergedTrackSource_,mergedTrajHandle);
00061   edm::Handle< TrajTrackAssociationCollection >  mergedTrajTrackHandle;
00062   iEvent.getByLabel(mergedTrackSource_,mergedTrajTrackHandle);
00063 
00064   edm::Handle< std::vector<Trajectory> >  originalTrajHandle;
00065   iEvent.getByLabel(originalTrackSource_,originalTrajHandle);
00066   edm::Handle< TrajTrackAssociationCollection >  originalTrajTrackHandle;
00067   iEvent.getByLabel(originalTrackSource_,originalTrajTrackHandle);
00068 
00069   edm::Handle<edm::View<DuplicateRecord> > candidateHandle;
00070   iEvent.getByLabel(candidateSource_,candidateHandle);
00071 
00072   std::auto_ptr<std::vector<reco::Track> > out_generalTracks(new std::vector<reco::Track>());
00073   out_generalTracks->reserve(originalHandle->size());
00074   reco::TrackRefProd refTrks = iEvent.getRefBeforePut<reco::TrackCollection>();
00075   std::auto_ptr< std::vector<Trajectory> > outputTrajs = std::auto_ptr< std::vector<Trajectory> >(new std::vector<Trajectory>());
00076   outputTrajs->reserve(originalTrajHandle->size()+mergedTrajHandle->size());
00077   edm::RefProd< std::vector<Trajectory> > refTrajs;
00078   std::auto_ptr< TrajTrackAssociationCollection >  outputTTAss = std::auto_ptr< TrajTrackAssociationCollection >(new TrajTrackAssociationCollection());
00079   //std::auto_ptr< TrajectorySeedCollection > outputSeeds
00080 
00081   std::auto_ptr<reco::TrackExtraCollection> outputTrkExtras;
00082   reco::TrackExtraRefProd refTrkExtras;
00083   std::auto_ptr<TrackingRecHitCollection> outputTrkHits;
00084   TrackingRecHitRefProd refTrkHits;
00085   std::auto_ptr<TrajectorySeedCollection> outputSeeds;
00086   edm::RefProd< TrajectorySeedCollection > refTrajSeeds;
00087 
00088   const int rSize = (int)originalHandle->size();
00089   edm::RefToBase<TrajectorySeed> seedsRefs[rSize];
00090 
00091   if(copyExtras_){
00092     outputTrkExtras = std::auto_ptr<reco::TrackExtraCollection>(new reco::TrackExtraCollection);
00093     outputTrkExtras->reserve(originalHandle->size());
00094     refTrkExtras    = iEvent.getRefBeforePut<reco::TrackExtraCollection>();
00095     outputTrkHits   = std::auto_ptr<TrackingRecHitCollection>(new TrackingRecHitCollection);
00096     outputTrkHits->reserve(originalHandle->size()*25);
00097     refTrkHits      = iEvent.getRefBeforePut<TrackingRecHitCollection>();
00098     if (makeReKeyedSeeds_){
00099       outputSeeds = std::auto_ptr<TrajectorySeedCollection>(new TrajectorySeedCollection);
00100       outputSeeds->reserve(originalHandle->size());
00101       refTrajSeeds = iEvent.getRefBeforePut<TrajectorySeedCollection>();
00102     }
00103   }
00104 
00105   //match new tracks to their candidates
00106   std::vector<std::pair<int,DuplicateRecord> > matches;
00107   for(int i = 0; i < (int)candidateHandle->size(); i++){
00108     DuplicateRecord duplicateRecord = candidateHandle->at(i);
00109     int matchTrack = matchCandidateToTrack(duplicateRecord.first,mergedHandle);
00110     if(matchTrack < 0)continue;
00111     if( ChiSquaredProbability(mergedHandle->at(matchTrack).chi2(),mergedHandle->at(matchTrack).ndof()) < minTrkProbCut_)continue;
00112     unsigned int dHits = (duplicateRecord.first.recHits().second - duplicateRecord.first.recHits().first) - mergedHandle->at(matchTrack).recHitsSize();
00113     if(dHits > diffHitsCut_)continue;
00114     matches.push_back(std::pair<int,DuplicateRecord>(matchTrack,duplicateRecord));
00115   }
00116 
00117   //check for candidates/tracks that share merged tracks, select minimum chi2, remove the rest
00118   std::vector<std::pair<int,DuplicateRecord> >::iterator matchIter0 = matches.begin();
00119   std::vector<std::pair<int,DuplicateRecord> >::iterator matchIter1;
00120   while(matchIter0 != matches.end()){
00121     double nchi2 = (mergedHandle->at((*matchIter0).first)).normalizedChi2();
00122     bool advance = true;
00123     for(matchIter1 = matchIter0+1; matchIter1 != matches.end(); matchIter1++){
00124       if((*matchIter1).second.second.first.seedRef() == (*matchIter0).second.second.first.seedRef() || (*matchIter1).second.second.first.seedRef() == (*matchIter0).second.second.second.seedRef() || (*matchIter1).second.second.second.seedRef() == (*matchIter0).second.second.first.seedRef() || (*matchIter1).second.second.second.seedRef() == (*matchIter0).second.second.second.seedRef()){
00125         double nchi2_1 = (mergedHandle->at((*matchIter1).first)).normalizedChi2();
00126         advance = false;
00127         if(nchi2_1 < nchi2){
00128           matches.erase(matchIter0);
00129         }else{
00130           matches.erase(matchIter1);
00131         }
00132         break;
00133       }
00134     }
00135     if(advance)matchIter0++;
00136   }
00137 
00138   //add the good merged tracks to the output list, remove input tracks
00139   std::vector<reco::Track> inputTracks;
00140 
00141   refTrajs = iEvent.getRefBeforePut< std::vector<Trajectory> >();
00142 
00143   for(matchIter0 = matches.begin(); matchIter0 != matches.end(); matchIter0++){
00144     reco::Track inTrk1 = (*matchIter0).second.second.first;
00145     reco::Track inTrk2 = (*matchIter0).second.second.second;
00146     reco::TrackBase::TrackAlgorithm newTrkAlgo = std::min(inTrk1.algo(),inTrk2.algo());
00147     int combinedQualityMask = (inTrk1.qualityMask() | inTrk2.qualityMask());
00148     inputTracks.push_back(inTrk1);
00149     inputTracks.push_back(inTrk2);
00150     out_generalTracks->push_back(mergedHandle->at((*matchIter0).first));
00151     reco::TrackRef curTrackRef = reco::TrackRef(refTrks, out_generalTracks->size() - 1);
00152     out_generalTracks->back().setAlgorithm(newTrkAlgo);
00153     out_generalTracks->back().setQualityMask(combinedQualityMask);
00154     out_generalTracks->back().setQuality(qualityToSet_);
00155     edm::RefToBase<TrajectorySeed> origSeedRef;
00156     if(copyExtras_){
00157       const reco::Track track = mergedHandle->at((*matchIter0).first);
00158       origSeedRef = track.seedRef();
00159       //creating a seed with rekeyed clusters if required
00160       if (makeReKeyedSeeds_){
00161         bool doRekeyOnThisSeed=false;
00162         
00163         edm::InputTag clusterRemovalInfos("");
00164         //grab on of the hits of the seed
00165         if (origSeedRef->nHits()!=0){
00166           TrajectorySeed::const_iterator firstHit=origSeedRef->recHits().first;
00167           const TrackingRecHit *hit = &*firstHit;
00168           if (firstHit->isValid()){
00169             edm::ProductID  pID=clusterProductB(hit);
00170             // the cluster collection either produced a removalInfo or mot
00171             //get the clusterremoval info from the provenance: will rekey if this is found
00172             edm::Handle<reco::ClusterRemovalInfo> CRIh;
00173             edm::Provenance prov=iEvent.getProvenance(pID);
00174             clusterRemovalInfos=edm::InputTag(prov.moduleLabel(),
00175                                               prov.productInstanceName(),
00176                                               prov.processName());
00177             doRekeyOnThisSeed=iEvent.getByLabel(clusterRemovalInfos,CRIh);
00178           }//valid hit
00179         }//nhit!=0
00180         
00181         if (doRekeyOnThisSeed && !(clusterRemovalInfos==edm::InputTag(""))) {
00182           ClusterRemovalRefSetter refSetter(iEvent,clusterRemovalInfos);
00183           TrajectorySeed::recHitContainer  newRecHitContainer;
00184           newRecHitContainer.reserve(origSeedRef->nHits());
00185           TrajectorySeed::const_iterator iH=origSeedRef->recHits().first;
00186           TrajectorySeed::const_iterator iH_end=origSeedRef->recHits().second;
00187           for (;iH!=iH_end;++iH){
00188             newRecHitContainer.push_back(*iH);
00189             refSetter.reKey(&newRecHitContainer.back());
00190           }
00191           outputSeeds->push_back( TrajectorySeed( origSeedRef->startingState(),
00192                                                   newRecHitContainer,
00193                                                   origSeedRef->direction()));
00194         }
00195         //doRekeyOnThisSeed=true
00196         else{
00197           //just copy the one we had before
00198           outputSeeds->push_back( TrajectorySeed(*origSeedRef));
00199           }
00200         edm::Ref<TrajectorySeedCollection> pureRef(refTrajSeeds, outputSeeds->size()-1);
00201         origSeedRef=edm::RefToBase<TrajectorySeed>( pureRef);
00202       }//creating a new seed and rekeying it rechit clusters.
00203       // Fill TrackExtra collection
00204       outputTrkExtras->push_back( reco::TrackExtra( 
00205                                                    track.outerPosition(), track.outerMomentum(), track.outerOk(),
00206                                                    track.innerPosition(), track.innerMomentum(), track.innerOk(),
00207                                                    track.outerStateCovariance(), track.outerDetId(),
00208                                                    track.innerStateCovariance(), track.innerDetId(),
00209                                                    track.seedDirection(), origSeedRef ) );
00210       seedsRefs[(*matchIter0).first]=origSeedRef;
00211       out_generalTracks->back().setExtra( reco::TrackExtraRef( refTrkExtras, outputTrkExtras->size() - 1) );
00212       reco::TrackExtra & tx = outputTrkExtras->back();
00213       tx.setResiduals(track.residuals());
00214       // fill TrackingRecHits
00215       unsigned nh1=track.recHitsSize();
00216       for ( unsigned ih=0; ih<nh1; ++ih ) { 
00217           //const TrackingRecHit*hit=&((*(track->recHit(ih))));
00218         outputTrkHits->push_back( track.recHit(ih)->clone() );
00219         tx.add( TrackingRecHitRef( refTrkHits, outputTrkHits->size() - 1) );
00220       }
00221     }
00222     edm::Ref< std::vector<Trajectory> > trajRef(mergedTrajHandle, (*matchIter0).first);
00223     TrajTrackAssociationCollection::const_iterator match = mergedTrajTrackHandle->find(trajRef);
00224     if (match != mergedTrajTrackHandle->end()) {
00225         if(curTrackRef.isNonnull()){
00226           outputTrajs->push_back( *trajRef );
00227           if (copyExtras_ && makeReKeyedSeeds_)
00228             outputTrajs->back().setSeedRef( origSeedRef );
00229           outputTTAss->insert(edm::Ref< std::vector<Trajectory> >(refTrajs, outputTrajs->size() - 1),curTrackRef );
00230         }
00231     }
00232   }
00233 
00234   for(int i = 0; i < (int)originalHandle->size(); i++){
00235     bool good = true;
00236     reco::Track origTrack = originalHandle->at(i);
00237     for(int j = 0; j < (int)inputTracks.size() && good; j++){
00238       reco::Track inputTrack = inputTracks[j];
00239       if(origTrack.seedRef() != inputTrack.seedRef())continue;
00240       if(origTrack.charge() != inputTrack.charge())continue;
00241       if(origTrack.momentum() != inputTrack.momentum())continue;
00242       if(origTrack.referencePoint() != inputTrack.referencePoint())continue;
00243       good = false;
00244     }
00245 
00246     if(good){
00247       out_generalTracks->push_back(origTrack);
00248       reco::TrackRef curTrackRef = reco::TrackRef(refTrks, out_generalTracks->size() - 1);
00249       edm::RefToBase<TrajectorySeed> origSeedRef;
00250 
00251       if(copyExtras_){
00252         const reco::Track track = origTrack;
00253         origSeedRef = track.seedRef();
00254         //creating a seed with rekeyed clusters if required
00255         if (makeReKeyedSeeds_){
00256           bool doRekeyOnThisSeed=false;
00257           
00258           edm::InputTag clusterRemovalInfos("");
00259           //grab on of the hits of the seed
00260           if (origSeedRef->nHits()!=0){
00261             TrajectorySeed::const_iterator firstHit=origSeedRef->recHits().first;
00262             const TrackingRecHit *hit = &*firstHit;
00263             if (firstHit->isValid()){
00264               edm::ProductID  pID=clusterProductB(hit);
00265               // the cluster collection either produced a removalInfo or mot
00266               //get the clusterremoval info from the provenance: will rekey if this is found
00267               edm::Handle<reco::ClusterRemovalInfo> CRIh;
00268               edm::Provenance prov=iEvent.getProvenance(pID);
00269               clusterRemovalInfos=edm::InputTag(prov.moduleLabel(),
00270                                                 prov.productInstanceName(),
00271                                                 prov.processName());
00272               doRekeyOnThisSeed=iEvent.getByLabel(clusterRemovalInfos,CRIh);
00273             }//valid hit
00274           }//nhit!=0
00275           
00276           if (doRekeyOnThisSeed && !(clusterRemovalInfos==edm::InputTag(""))) {
00277             ClusterRemovalRefSetter refSetter(iEvent,clusterRemovalInfos);
00278             TrajectorySeed::recHitContainer  newRecHitContainer;
00279             newRecHitContainer.reserve(origSeedRef->nHits());
00280             TrajectorySeed::const_iterator iH=origSeedRef->recHits().first;
00281             TrajectorySeed::const_iterator iH_end=origSeedRef->recHits().second;
00282             for (;iH!=iH_end;++iH){
00283               newRecHitContainer.push_back(*iH);
00284               refSetter.reKey(&newRecHitContainer.back());
00285             }
00286             outputSeeds->push_back( TrajectorySeed( origSeedRef->startingState(),
00287                                                     newRecHitContainer,
00288                                                     origSeedRef->direction()));
00289           }
00290           //doRekeyOnThisSeed=true
00291           else{
00292             //just copy the one we had before
00293             outputSeeds->push_back( TrajectorySeed(*origSeedRef));
00294           }
00295           edm::Ref<TrajectorySeedCollection> pureRef(refTrajSeeds, outputSeeds->size()-1);
00296           origSeedRef=edm::RefToBase<TrajectorySeed>( pureRef);
00297         }//creating a new seed and rekeying it rechit clusters.
00298         
00299         // Fill TrackExtra collection
00300         outputTrkExtras->push_back( reco::TrackExtra( 
00301                                                      track.outerPosition(), track.outerMomentum(), track.outerOk(),
00302                                                      track.innerPosition(), track.innerMomentum(), track.innerOk(),
00303                                                      track.outerStateCovariance(), track.outerDetId(),
00304                                                      track.innerStateCovariance(), track.innerDetId(),
00305                                                      track.seedDirection(), origSeedRef ) );
00306         seedsRefs[i]=origSeedRef;
00307         out_generalTracks->back().setExtra( reco::TrackExtraRef( refTrkExtras, outputTrkExtras->size() - 1) );
00308         reco::TrackExtra & tx = outputTrkExtras->back();
00309         tx.setResiduals(track.residuals());
00310         
00311         // fill TrackingRecHits
00312         unsigned nh1=track.recHitsSize();
00313         for ( unsigned ih=0; ih<nh1; ++ih ) { 
00314           //const TrackingRecHit*hit=&((*(track->recHit(ih))));
00315           outputTrkHits->push_back( track.recHit(ih)->clone() );
00316           tx.add( TrackingRecHitRef( refTrkHits, outputTrkHits->size() - 1) );
00317         }
00318         
00319       }
00320 
00321       edm::Ref< std::vector<Trajectory> > trajRef(originalTrajHandle, i);
00322       TrajTrackAssociationCollection::const_iterator match = originalTrajTrackHandle->find(trajRef);
00323       if (match != originalTrajTrackHandle->end()) {
00324         if(curTrackRef.isNonnull()){
00325           outputTrajs->push_back( *trajRef );
00326           if (copyExtras_ && makeReKeyedSeeds_)
00327             outputTrajs->back().setSeedRef( origSeedRef );
00328           outputTTAss->insert(edm::Ref< std::vector<Trajectory> >(refTrajs, outputTrajs->size() - 1),curTrackRef );
00329         }
00330       }
00331     }
00332   }
00333 
00334   iEvent.put(out_generalTracks);
00335   if (copyExtras_) {
00336     iEvent.put(outputTrkExtras);
00337     iEvent.put(outputTrkHits);
00338     if (makeReKeyedSeeds_)
00339       iEvent.put(outputSeeds);
00340   }
00341   iEvent.put(outputTrajs);
00342   iEvent.put(outputTTAss);
00343 }
00344 
00345 //------------------------------------------------------------
00346 //------------------------------------------------------------
00347 //------------------------------------------------------------
00348 //------------------------------------------------------------
00349 int DuplicateListMerger::matchCandidateToTrack(TrackCandidate candidate, edm::Handle<edm::View<reco::Track> > tracks){
00350   int track = -1;
00351   std::vector<int> rawIds;
00352   RecHitContainer::const_iterator candBegin = candidate.recHits().first;
00353   RecHitContainer::const_iterator candEnd = candidate.recHits().second;
00354   for(; candBegin != candEnd; candBegin++){
00355     rawIds.push_back((*(candBegin)).rawId());
00356   }
00357  
00358 
00359   for(int i = 0; i < (int)tracks->size() && track < 0;i++){
00360     if((tracks->at(i)).seedRef() != candidate.seedRef())continue;
00361     int match = 0;
00362     trackingRecHit_iterator trackRecBegin = tracks->at(i).recHitsBegin();
00363     trackingRecHit_iterator trackRecEnd = tracks->at(i).recHitsEnd();
00364     for(;trackRecBegin != trackRecEnd; trackRecBegin++){
00365       if(std::find(rawIds.begin(),rawIds.end(),(*(trackRecBegin)).get()->rawId()) != rawIds.end())match++;
00366     }
00367     if(match != (int)tracks->at(i).recHitsSize())continue;
00368     track = i;
00369   }
00370 
00371   return track;
00372 }
00373 //------------------------------------------------------------
00374 //------------------------------------------------------------
00375 //------------------------------------------------------------
00376 //------------------------------------------------------------