CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/RecoMuon/TrackingTools/src/MuonTrajectoryCleaner.cc

Go to the documentation of this file.
00001 
00008 #include "RecoMuon/TrackingTools/interface/MuonTrajectoryCleaner.h"
00009 #include "DataFormats/TrackReco/interface/Track.h"
00010 #include "DataFormats/Math/interface/deltaR.h"
00011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00012 #include "DataFormats/Common/interface/RefToBase.h"
00013 #include "DataFormats/Common/interface/AssociationMap.h"
00014 #include "DataFormats/MuonSeed/interface/L2MuonTrajectorySeedCollection.h"
00015 
00016 using namespace std; 
00017 
00018 typedef edm::AssociationMap<edm::OneToMany<L2MuonTrajectorySeedCollection, L2MuonTrajectorySeedCollection> > L2SeedAssoc;
00019 
00020 void MuonTrajectoryCleaner::clean(TrajectoryContainer& trajC, edm::Event& event) {
00021   const std::string metname = "Muon|RecoMuon|MuonTrajectoryCleaner";
00022 
00023   LogTrace(metname) << "Muon Trajectory Cleaner called" << endl;
00024 
00025   TrajectoryContainer::iterator iter, jter;
00026   Trajectory::DataContainer::const_iterator m1, m2;
00027 
00028   if ( trajC.size() < 1 ) return;
00029   
00030   LogTrace(metname) << "Number of trajectories in the container: " <<trajC.size()<< endl;
00031 
00032   int i(0), j(0);
00033   int match(0);
00034 
00035   // Map between chosen seed and ghost seeds (for trigger)
00036   map<int, vector<int> > seedmap;
00037 
00038   // CAVEAT: vector<bool> is not a vector, its elements are not addressable!
00039   // This is fine as long as only operator [] is used as in this case.
00040   // cf. par 16.3.11
00041   vector<bool> mask(trajC.size(),true);
00042   
00043   TrajectoryContainer result;
00044   
00045   for ( iter = trajC.begin(); iter != trajC.end(); iter++ ) {
00046     if ( !mask[i] ) { i++; continue; }
00047     if ( !(*iter)->isValid() || (*iter)->empty() ) {mask[i] = false; i++; continue; }
00048     if(seedmap.count(i)==0) seedmap[i].push_back(i);
00049     const Trajectory::DataContainer& meas1 = (*iter)->measurements();
00050     j = i+1;
00051     bool skipnext=false;
00052     for ( jter = iter+1; jter != trajC.end(); jter++ ) {
00053       if ( !mask[j] ) { j++; continue; }
00054       if(seedmap.count(j)==0) seedmap[j].push_back(j);
00055       const Trajectory::DataContainer& meas2 = (*jter)->measurements();
00056       match = 0;
00057       for ( m1 = meas1.begin(); m1 != meas1.end(); m1++ ) {
00058         for ( m2 = meas2.begin(); m2 != meas2.end(); m2++ ) {
00059           if ( ( (*m1).recHit()->globalPosition() - (*m2).recHit()->globalPosition()).mag()< 10e-5 ) match++;
00060         }
00061       }
00062       
00063 
00064       // FIXME Set Boff/on via cfg!
00065       double chi2_dof_i = (*iter)->ndof() > 0 ? (*iter)->chiSquared()/(*iter)->ndof() : (*iter)->chiSquared()/1e-10;
00066       double chi2_dof_j = (*jter)->ndof() > 0 ? (*jter)->chiSquared()/(*jter)->ndof() : (*jter)->chiSquared()/1e-10;
00067 
00068       LogTrace(metname) 
00069         << " MuonTrajectoryCleaner: trajC " 
00070         << i << "(pT="<<(*iter)->lastMeasurement().updatedState().globalMomentum().perp() << " GeV)"
00071         << " chi2/nDOF = " << (*iter)->chiSquared() << "/" << (*iter)->ndof() 
00072         << " (RH=" << (*iter)->foundHits() << ") = " << chi2_dof_i
00073         << " vs trajC " 
00074         << j << "(pT="<<(*jter)->lastMeasurement().updatedState().globalMomentum().perp() << " GeV)"
00075         << " chi2/nRH = " << (*jter)->chiSquared() << "/" <<  (*jter)->ndof() 
00076         << " (RH=" << (*jter)->foundHits() << ") = " << chi2_dof_j
00077         << " Shared RecHits: " << match; 
00078 
00079       int hit_diff =  (*iter)->foundHits() - (*jter)->foundHits() ;       
00080       // If there are matches, reject the worst track
00081       if ( match > 0 ) {
00082         // If the difference of # of rechits is less than 4, compare the chi2/ndf
00083         if ( abs(hit_diff) <= 4  ) {
00084 
00085           double minPt = 3.5;
00086           double dPt = 7.;  // i.e. considering 10% (conservative!) resolution at minPt it is ~ 10 sigma away from the central value
00087 
00088           double maxFraction = 0.95;
00089 
00090           double fraction = (2.*match)/((*iter)->measurements().size()+(*jter)->measurements().size());
00091           int belowLimit = 0;
00092           int above = 0;
00093 
00094           if((*jter)->lastMeasurement().updatedState().globalMomentum().perp() <= minPt) ++belowLimit; 
00095           if((*iter)->lastMeasurement().updatedState().globalMomentum().perp() <= minPt) ++belowLimit; 
00096          
00097           if((*jter)->lastMeasurement().updatedState().globalMomentum().perp() >= dPt) ++above; 
00098           if((*iter)->lastMeasurement().updatedState().globalMomentum().perp() >= dPt) ++above; 
00099           
00100 
00101           if(fraction >= maxFraction && belowLimit == 1 && above == 1){
00102             if((*iter)->lastMeasurement().updatedState().globalMomentum().perp() < minPt){
00103               mask[i] = false;
00104               skipnext=true;
00105               seedmap[j].insert(seedmap[j].end(), seedmap[i].begin(), seedmap[i].end());
00106               seedmap.erase(i);
00107               LogTrace(metname) << "Trajectory # " << i << " (pT="<<(*iter)->lastMeasurement().updatedState().globalMomentum().perp() 
00108                                 << " GeV) rejected because it has too low pt";
00109             } 
00110             else {
00111               mask[j] = false;
00112               seedmap[i].insert(seedmap[i].end(), seedmap[j].begin(), seedmap[j].end());
00113               seedmap.erase(j);
00114               LogTrace(metname) << "Trajectory # " << j << " (pT="<<(*jter)->lastMeasurement().updatedState().globalMomentum().perp() 
00115                                 << " GeV) rejected because it has too low pt";
00116             }
00117           }
00118           else{   
00119             if (chi2_dof_i  > chi2_dof_j) {
00120               mask[i] = false;
00121               skipnext=true;
00122               seedmap[j].insert(seedmap[j].end(), seedmap[i].begin(), seedmap[i].end());
00123               seedmap.erase(i);
00124               LogTrace(metname) << "Trajectory # " << i << " (pT="<<(*iter)->lastMeasurement().updatedState().globalMomentum().perp() << " GeV) rejected";
00125             }
00126             else{
00127               mask[j] = false;
00128               seedmap[i].insert(seedmap[i].end(), seedmap[j].begin(), seedmap[j].end());
00129               seedmap.erase(j);
00130               LogTrace(metname) << "Trajectory # " << j << " (pT="<<(*jter)->lastMeasurement().updatedState().globalMomentum().perp() << " GeV) rejected";
00131             }
00132           }
00133         }
00134         else { // different number of hits
00135           if ( hit_diff < 0 ) {
00136             mask[i] = false;
00137             skipnext=true;
00138             seedmap[j].insert(seedmap[j].end(), seedmap[i].begin(), seedmap[i].end());
00139             seedmap.erase(i);
00140             LogTrace(metname) << "Trajectory # " << i << " (pT="<<(*iter)->lastMeasurement().updatedState().globalMomentum().perp() << " GeV) rejected";
00141           }
00142           else { 
00143             mask[j] = false;
00144             seedmap[i].insert(seedmap[i].end(), seedmap[j].begin(), seedmap[j].end());
00145             seedmap.erase(j);
00146             LogTrace(metname) << "Trajectory # " << j << " (pT="<<(*jter)->lastMeasurement().updatedState().globalMomentum().perp() << " GeV) rejected";
00147           }
00148         } 
00149       }
00150       if(skipnext) break;
00151       j++;
00152     }
00153     i++;
00154     if(skipnext) continue;
00155   }
00156   
00157 
00158   // Association map between the seed of the chosen trajectory and the seeds of the discarded trajectories
00159   if(reportGhosts_) {
00160     LogTrace(metname) << " Creating map between chosen seed and ghost seeds." << std::endl;
00161 
00162     auto_ptr<L2SeedAssoc> seedToSeedsMap(new L2SeedAssoc);
00163     int seedcnt(0);
00164 
00165     for(map<int, vector<int> >::iterator itmap=seedmap.begin(); itmap!=seedmap.end(); ++itmap, ++seedcnt) {
00166       edm::RefToBase<TrajectorySeed> tmpSeedRef1 = trajC[(*itmap).first]->seedRef();
00167       edm::Ref<L2MuonTrajectorySeedCollection> tmpL2SeedRef1 = tmpSeedRef1.castTo<edm::Ref<L2MuonTrajectorySeedCollection> >();
00168 
00169       int tmpsize=(*itmap).second.size();
00170       for(int cnt=0; cnt!=tmpsize; ++cnt) {
00171         edm::RefToBase<TrajectorySeed> tmpSeedRef2 = trajC[(*itmap).second[cnt]]->seedRef();
00172         edm::Ref<L2MuonTrajectorySeedCollection> tmpL2SeedRef2 = tmpSeedRef2.castTo<edm::Ref<L2MuonTrajectorySeedCollection> >();
00173         seedToSeedsMap->insert(tmpL2SeedRef1, tmpL2SeedRef2);
00174       }
00175     }
00176 
00177     event.put(seedToSeedsMap, "");
00178 
00179   } // end if(reportGhosts_)
00180 
00181   i = 0;
00182   for ( iter = trajC.begin(); iter != trajC.end(); iter++ ) {
00183     if ( mask[i] ){
00184       result.push_back(*iter);
00185       LogTrace(metname) << "Keep trajectory with pT = " << (*iter)->lastMeasurement().updatedState().globalMomentum().perp() << " GeV";
00186     }
00187     else delete *iter;
00188     i++;
00189   }
00190 
00191   trajC.clear();
00192   trajC = result;
00193 }
00194 
00195 //
00196 // clean CandidateContainer
00197 //
00198 void MuonTrajectoryCleaner::clean(CandidateContainer& candC){ 
00199   const std::string metname = "Muon|RecoMuon|MuonTrajectoryCleaner";
00200 
00201   LogTrace(metname) << "Muon Trajectory Cleaner called" << endl;
00202 
00203   if ( candC.size() < 2 ) return;
00204 
00205   CandidateContainer::iterator iter, jter;
00206   Trajectory::DataContainer::const_iterator m1, m2;
00207 
00208   const float deltaEta = 0.01;
00209   const float deltaPhi = 0.01;
00210   const float deltaPt  = 1.0;
00211   
00212   LogTrace(metname) << "Number of muon candidates in the container: " <<candC.size()<< endl;
00213 
00214   int i(0), j(0);
00215   int match(0);
00216   bool directionMatch = false;
00217 
00218   // CAVEAT: vector<bool> is not a vector, its elements are not addressable!
00219   // This is fine as long as only operator [] is used as in this case.
00220   // cf. par 16.3.11
00221   vector<bool> mask(candC.size(),true);
00222   
00223   CandidateContainer result;
00224   
00225   for ( iter = candC.begin(); iter != candC.end(); iter++ ) {
00226     if ( !mask[i] ) { i++; continue; }
00227     const Trajectory::DataContainer& meas1 = (*iter)->trajectory()->measurements();
00228     j = i+1;
00229     bool skipnext=false;
00230 
00231     TrajectoryStateOnSurface innerTSOS;
00232 
00233     if ((*iter)->trajectory()->direction() == alongMomentum) {
00234       innerTSOS = (*iter)->trajectory()->firstMeasurement().updatedState();
00235     } 
00236     else if ((*iter)->trajectory()->direction() == oppositeToMomentum) { 
00237       innerTSOS = (*iter)->trajectory()->lastMeasurement().updatedState();
00238          }
00239     if ( !(innerTSOS.isValid()) ) continue;
00240     float pt1 = innerTSOS.globalMomentum().perp();
00241     float eta1 = innerTSOS.globalMomentum().eta();
00242     float phi1 = innerTSOS.globalMomentum().phi();
00243 
00244     for ( jter = iter+1; jter != candC.end(); jter++ ) {
00245       if ( !mask[j] ) { j++; continue; }
00246       directionMatch = false;
00247       const Trajectory::DataContainer& meas2 = (*jter)->trajectory()->measurements();
00248       match = 0;
00249       for ( m1 = meas1.begin(); m1 != meas1.end(); m1++ ) {
00250         for ( m2 = meas2.begin(); m2 != meas2.end(); m2++ ) {
00251           if ( (*m1).recHit()->isValid() && (*m2).recHit()->isValid() ) 
00252             if ( ( (*m1).recHit()->globalPosition() - (*m2).recHit()->globalPosition()).mag()< 10e-5 ) match++;
00253         }
00254       }
00255       
00256       LogTrace(metname) 
00257         << " MuonTrajectoryCleaner: candC " << i << " chi2/nRH = " 
00258         << (*iter)->trajectory()->chiSquared() << "/" << (*iter)->trajectory()->foundHits() <<
00259         " vs trajC " << j << " chi2/nRH = " << (*jter)->trajectory()->chiSquared() <<
00260         "/" << (*jter)->trajectory()->foundHits() << " Shared RecHits: " << match;
00261 
00262       TrajectoryStateOnSurface innerTSOS2;       
00263       if ((*jter)->trajectory()->direction() == alongMomentum) {
00264         innerTSOS2 = (*jter)->trajectory()->firstMeasurement().updatedState();
00265       }
00266       else if ((*jter)->trajectory()->direction() == oppositeToMomentum) {
00267         innerTSOS2 = (*jter)->trajectory()->lastMeasurement().updatedState();
00268       }
00269       if ( !(innerTSOS2.isValid()) ) continue;
00270 
00271       float pt2 = innerTSOS2.globalMomentum().perp();
00272       float eta2 = innerTSOS2.globalMomentum().eta();
00273       float phi2 = innerTSOS2.globalMomentum().phi();
00274 
00275       float deta(fabs(eta1-eta2));
00276       float dphi(fabs(Geom::Phi<float>(phi1)-Geom::Phi<float>(phi2)));
00277       float dpt(fabs(pt1-pt2));
00278       if ( dpt < deltaPt && deta < deltaEta && dphi < deltaPhi ) {
00279         directionMatch = true;
00280         LogTrace(metname)
00281         << " MuonTrajectoryCleaner: candC " << i<<" and "<<j<< " direction matched: "
00282         <<innerTSOS.globalMomentum()<<" and " <<innerTSOS2.globalMomentum();
00283       }
00284 
00285       // If there are matches, reject the worst track
00286       bool hitsMatch = ((match > 0) && ((match/((*iter)->trajectory()->foundHits()) > 0.25) || (match/((*jter)->trajectory()->foundHits()) > 0.25))) ? true : false;
00287       bool tracksMatch = 
00288       ( ( (*iter)->trackerTrack() == (*jter)->trackerTrack() ) && 
00289         ( deltaR<double>((*iter)->muonTrack()->eta(),(*iter)->muonTrack()->phi(), (*jter)->muonTrack()->eta(),(*jter)->muonTrack()->phi()) < 0.2) );
00290 
00291       //if ( ( tracksMatch ) || (hitsMatch > 0) || directionMatch ) {
00292                         if ( ( tracksMatch ) || (hitsMatch > 0) ) { 
00293         if (  (*iter)->trajectory()->foundHits() == (*jter)->trajectory()->foundHits() ) {
00294           if ( (*iter)->trajectory()->chiSquared() > (*jter)->trajectory()->chiSquared() ) {
00295             mask[i] = false;
00296             skipnext=true;
00297           }
00298           else mask[j] = false;
00299         }
00300         else { // different number of hits
00301           if ( (*iter)->trajectory()->foundHits() < (*jter)->trajectory()->foundHits() ) {
00302             mask[i] = false;
00303             skipnext=true;
00304           }
00305           else mask[j] = false;
00306         }
00307       }
00308       if(skipnext) break;
00309       j++;
00310     }
00311     i++;
00312     if(skipnext) continue;
00313   }
00314   
00315   i = 0;
00316   for ( iter = candC.begin(); iter != candC.end(); iter++ ) {
00317     if ( mask[i] ) {
00318        result.push_back(*iter);
00319     } else {
00320        delete (*iter)->trajectory();
00321        delete (*iter)->trackerTrajectory();
00322        delete *iter;
00323     } 
00324     i++;
00325   }
00326   
00327   candC.clear();
00328   candC = result;
00329 }
00330