00001 // 00002 // Package: RecoTracker/RoadSearchTrackListCleaner 00003 // Class: RoadSearchHelixMaker 00004 // 00005 // Description: TrackListCleaner 00006 // 00007 // Original Author: Steve Wagner, stevew@pizero.colorado.edu 00008 // Created: Sat Jan 14 22:00:00 UTC 2006 00009 // 00010 // $Author: stevew $ 00011 // $Date: 2007/07/16 20:25:35 $ 00012 // $Revision: 1.6 $ 00013 // 00014 00015 #include <memory> 00016 #include <string> 00017 #include <iostream> 00018 #include <cmath> 00019 #include <vector> 00020 00021 #include "RecoTracker/RoadSearchHelixMaker/interface/RoadSearchTrackListCleaner.h" 00022 00023 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2DCollection.h" 00024 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2DCollection.h" 00025 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h" 00026 #include "DataFormats/TrackCandidate/interface/TrackCandidateCollection.h" 00027 00028 #include "DataFormats/Common/interface/Handle.h" 00029 #include "FWCore/Framework/interface/ESHandle.h" 00030 #include "FWCore/Framework/interface/EventSetup.h" 00031 00032 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h" 00033 #include "Geometry/CommonDetUnit/interface/GeomDet.h" 00034 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h" 00035 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h" 00036 00037 #include "DataFormats/TrackReco/interface/TrackFwd.h" 00038 #include "DataFormats/TrackReco/interface/Track.h" 00039 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h" 00040 #include "DataFormats/TrackReco/interface/TrackExtra.h" 00041 #include "DataFormats/TrackReco/src/classes.h" 00042 00043 namespace cms 00044 { 00045 00046 RoadSearchTrackListCleaner::RoadSearchTrackListCleaner(edm::ParameterSet const& conf) : 00047 conf_(conf) 00048 { 00049 produces<reco::TrackCollection>(); 00050 // produces<reco::TrackExtraCollection>(); 00051 } 00052 00053 00054 // Virtual destructor needed. 00055 RoadSearchTrackListCleaner::~RoadSearchTrackListCleaner() { } 00056 00057 // Functions that gets called by framework every event 00058 void RoadSearchTrackListCleaner::produce(edm::Event& e, const edm::EventSetup& es) 00059 { 00060 // retrieve producer name of input SiStripRecHit2DCollection 00061 std::string trackProducer = conf_.getParameter<std::string>("TrackProducer"); 00062 double maxNormalizedChisq = conf_.getParameter<double>("MaxNormalizedChisq"); 00063 double minPT = conf_.getParameter<double>("MinPT"); 00064 unsigned int minFound = (unsigned int)conf_.getParameter<int>("MinFound"); 00065 00066 // 00067 // extract tracker geometry 00068 // 00069 edm::ESHandle<TrackerGeometry> theG; 00070 es.get<TrackerDigiGeometryRecord>().get(theG); 00071 00072 // using namespace reco; 00073 00074 // get Inputs 00075 edm::Handle<reco::TrackCollection> trackCollection; 00076 e.getByLabel(trackProducer, trackCollection); 00077 00078 const reco::TrackCollection tC = *(trackCollection.product()); 00079 00080 //std::cout << "Reconstructed "<< tC.size() << " tracks" << std::endl ; 00081 00082 // Step B: create empty output collection 00083 std::auto_ptr<reco::TrackCollection> output(new reco::TrackCollection); 00084 00085 // 00086 // no input tracks 00087 // 00088 00089 if ( tC.empty() ){ 00090 // LogDebug("RoadSearch") << "Found " << output.size() << " clouds."; 00091 e.put(output); 00092 return; 00093 } 00094 00095 // 00096 // quality cuts first 00097 // 00098 std::vector<int> selected; for (unsigned int i=0; i<tC.size(); ++i){selected.push_back(1);} 00099 00100 int i=-1; 00101 for (reco::TrackCollection::const_iterator track=tC.begin(); track!=tC.end(); track++){ 00102 i++; 00103 if ((short unsigned)track->ndof() < 1){ 00104 selected[i]=0; 00105 //std::cout << "Track "<< i << " rejected in TrackListCleaner; ndof() < 1" << std::endl ; 00106 continue; 00107 } 00108 if (track->normalizedChi2() > maxNormalizedChisq){ 00109 selected[i]=0; 00110 //std::cout << "Track "<< i << " rejected in TrackListCleaner; normalizedChi2() > maxNormalizedChisq " << track->normalizedChi2() << " " << maxNormalizedChisq << std::endl ; 00111 continue; 00112 } 00113 if (track->found() < minFound){ 00114 selected[i]=0; 00115 //std::cout << "Track "<< i << " rejected in TrackListCleaner; found() < minFound " << track->found() << " " << minFound << std::endl ; 00116 continue; 00117 } 00118 if (track->pt() < minPT){ 00119 selected[i]=0; 00120 //std::cout << "Track "<< i << " rejected in TrackListCleaner; pt() < minPT " << track->pt() << " " << minPT << std::endl ; 00121 continue; 00122 } 00123 }//end loop over tracks 00124 00125 // 00126 // > 1 track - try merging 00127 // 00128 if ( 1<tC.size() ){ 00129 int i=-1; 00130 for (reco::TrackCollection::const_iterator track=tC.begin(); track!=tC.end(); track++){ 00131 i++; 00132 //std::cout << "Track number "<< i << std::endl ; 00133 if (!selected[i])continue; 00134 int j=-1; 00135 for (reco::TrackCollection::const_iterator track2=tC.begin(); track2!=tC.end(); track2++){ 00136 j++; 00137 if ((j<=i)||(!selected[j])||(!selected[i]))continue; 00138 int noverlap=0; 00139 for (trackingRecHit_iterator it = track->recHitsBegin(); it != track->recHitsEnd(); it++){ 00140 if ((*it)->isValid()){ 00141 for (trackingRecHit_iterator jt = track2->recHitsBegin(); jt != track2->recHitsEnd(); jt++){ 00142 if ((*jt)->isValid()){ 00143 if (((*it)->geographicalId()==(*jt)->geographicalId())&&((*it)->localPosition().x()==(*jt)->localPosition().x()))noverlap++; 00144 } 00145 } 00146 } 00147 } 00148 float fi=float(noverlap)/float(track->recHitsSize()); float fj=float(noverlap)/float(track2->recHitsSize()); 00149 //std::cout << " trk1 trk2 nhits1 nhits2 nover " << i << " " << j << " " << track->recHitsSize() << " " << track2->recHitsSize() << " " << noverlap << " " << fi << " " << fj <<std::endl; 00150 if ((fi>0.66)||(fj>0.66)){ 00151 if (fi<fj){ 00152 selected[j]=0; 00153 //std::cout << " removing 2nd trk in pair " << std::endl; 00154 }else{ 00155 if (fi>fj){ 00156 selected[i]=0; 00157 //std::cout << " removing 1st trk in pair " << std::endl; 00158 }else{ 00159 //std::cout << " removing worst chisq in pair " << std::endl; 00160 if (track->chi2() > track2->chi2()){selected[i]=0;}else{selected[j]=0;} 00161 }//end fi > or = fj 00162 }//end fi < fj 00163 }//end got a duplicate 00164 }//end track2 loop 00165 }//end track loop 00166 }//end more than 1 track 00167 00168 // 00169 // output selected tracks - if any 00170 // 00171 i=-1; 00172 for (reco::TrackCollection::const_iterator track=tC.begin(); track!=tC.end(); track++){ 00173 i++; if (!selected[i])continue; 00174 reco::Track * theTrack = new reco::Track(track->chi2(), 00175 (short unsigned)track->ndof(), 00176 track->innerPosition(), 00177 track->innerMomentum(), 00178 track->charge(), 00179 track->innerStateCovariance()); 00180 //fill the TrackCollection 00181 reco::TrackExtraRef theTrackExtraRef=track->extra(); 00182 theTrack->setExtra(theTrackExtraRef); 00183 theTrack->setHitPattern((*theTrackExtraRef).recHits()); 00184 output->push_back(*theTrack); 00185 delete theTrack; 00186 }//end faux loop over tracks 00187 e.put(output); 00188 return; 00189 00190 }//end produce 00191 00192 }