CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/RecoEgamma/EgammaPhotonProducers/src/ConversionTrackMerger.cc

Go to the documentation of this file.
00001 //
00002 // Package:         RecoTracker/FinalTrackSelectors
00003 // Class:           ConversionTrackMerger
00004 // 
00005 // Description:     Merger for ConversionTracks, adapted from SimpleTrackListMerger
00006 //
00007 // Original Author: J.Bendavid
00008 //
00009 // $Author: bendavid $
00010 // $Date: 2010/11/22 02:02:08 $
00011 // $Revision: 1.2 $
00012 //
00013 
00014 #include <memory>
00015 #include <string>
00016 #include <iostream>
00017 #include <cmath>
00018 #include <vector>
00019 
00020 #include "RecoEgamma/EgammaPhotonProducers/interface/ConversionTrackMerger.h"
00021 
00022 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2DCollection.h"
00023 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2DCollection.h"
00024 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit1DCollection.h"
00025 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
00026 #include "DataFormats/TrackCandidate/interface/TrackCandidateCollection.h"
00027 
00028 #include "FWCore/Framework/interface/ESHandle.h"
00029 
00030 #include "FWCore/MessageLogger/interface/MessageLogger.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/src/classes.h"
00038 
00039 #include "RecoTracker/TrackProducer/interface/ClusterRemovalRefSetter.h"
00040     
00041   ConversionTrackMerger::ConversionTrackMerger(edm::ParameterSet const& conf) : 
00042     conf_(conf)
00043   {
00044     produces<reco::ConversionTrackCollection>();
00045    
00046   }
00047 
00048 
00049   // Virtual destructor needed.
00050   ConversionTrackMerger::~ConversionTrackMerger() { }  
00051 
00052   // Functions that gets called by framework every event
00053   void ConversionTrackMerger::produce(edm::Event& e, const edm::EventSetup& es)
00054   {
00055     // retrieve producer name of input TrackCollection(s)
00056     std::string trackProducer1 = conf_.getParameter<std::string>("TrackProducer1");
00057     std::string trackProducer2 = conf_.getParameter<std::string>("TrackProducer2");
00058 
00059     double shareFrac =  conf_.getParameter<double>("ShareFrac");
00060     bool allowFirstHitShare = conf_.getParameter<bool>("allowFirstHitShare");
00061     
00062     int outputPreferCollection = conf_.getParameter<int>("outputPreferCollection");
00063     int trackerOnlyPreferCollection = conf_.getParameter<int>("trackerOnlyPreferCollection");
00064     int arbitratedEcalSeededPreferCollection = conf_.getParameter<int>("arbitratedEcalSeededPreferCollection");    
00065     int arbitratedMergedPreferCollection = conf_.getParameter<int>("arbitratedMergedPreferCollection");
00066     int arbitratedMergedEcalGeneralPreferCollection = conf_.getParameter<int>("arbitratedMergedEcalGeneralPreferCollection");    
00067 
00068     // get Inputs 
00069     // if 1 input list doesn't exist, make an empty list, issue a warning, and continue
00070     // this allows ConversionTrackMerger to be used as a cleaner only if handed just one list
00071     // if both input lists don't exist, will issue 2 warnings and generate an empty output collection
00072     //
00073     const reco::ConversionTrackCollection *TC1 = 0;
00074     static const reco::ConversionTrackCollection s_empty1, s_empty2;
00075     edm::Handle<reco::ConversionTrackCollection> trackCollection1;
00076     e.getByLabel(trackProducer1, trackCollection1);
00077     if (trackCollection1.isValid()) {
00078       TC1 = trackCollection1.product();
00079       //std::cout << "1st collection " << trackProducer1 << " has "<< TC1->size() << " tracks" << std::endl ;
00080     } else {
00081       TC1 = &s_empty1;
00082       edm::LogWarning("ConversionTrackMerger") << "1st TrackCollection " << trackProducer1 << " not found; will only clean 2nd TrackCollection " << trackProducer2 ;
00083     }
00084     const reco::ConversionTrackCollection tC1 = *TC1;
00085 
00086     const reco::ConversionTrackCollection *TC2 = 0;
00087     edm::Handle<reco::ConversionTrackCollection> trackCollection2;
00088     e.getByLabel(trackProducer2, trackCollection2);
00089     if (trackCollection2.isValid()) {
00090       TC2 = trackCollection2.product();
00091       //std::cout << "2nd collection " << trackProducer2 << " has "<< TC2->size() << " tracks" << std::endl ;
00092     } else {
00093         TC2 = &s_empty2;
00094         edm::LogWarning("ConversionTrackMerger") << "2nd TrackCollection " << trackProducer2 << " not found; will only clean 1st TrackCollection " << trackProducer1 ;
00095     }
00096     const reco::ConversionTrackCollection tC2 = *TC2;
00097 
00098     // Step B: create empty output collection
00099     outputTrks = std::auto_ptr<reco::ConversionTrackCollection>(new reco::ConversionTrackCollection);
00100     int i;
00101 
00102     std::vector<int> selected1; for (unsigned int i=0; i<tC1.size(); ++i){selected1.push_back(1);}
00103     std::vector<int> selected2; for (unsigned int i=0; i<tC2.size(); ++i){selected2.push_back(1);}
00104 
00105    
00106    std::map<reco::ConversionTrackCollection::const_iterator, std::vector<const TrackingRecHit*> > rh1;
00107    std::map<reco::ConversionTrackCollection::const_iterator, std::vector<const TrackingRecHit*> > rh2;
00108    for (reco::ConversionTrackCollection::const_iterator track=tC1.begin(); track!=tC1.end(); ++track){
00109      trackingRecHit_iterator itB = track->track()->recHitsBegin();
00110      trackingRecHit_iterator itE = track->track()->recHitsEnd();
00111      for (trackingRecHit_iterator it = itB;  it != itE; ++it) { 
00112        const TrackingRecHit* hit = &(**it);
00113        rh1[track].push_back(hit);
00114      }
00115    }
00116    for (reco::ConversionTrackCollection::const_iterator track=tC2.begin(); track!=tC2.end(); ++track){
00117      trackingRecHit_iterator jtB = track->track()->recHitsBegin();
00118      trackingRecHit_iterator jtE = track->track()->recHitsEnd();
00119      for (trackingRecHit_iterator jt = jtB;  jt != jtE; ++jt) { 
00120        const TrackingRecHit* hit = &(**jt);
00121        rh2[track].push_back(hit);
00122      }
00123    }
00124 
00125    if ( (0<tC1.size())&&(0<tC2.size()) ){
00126     i=-1;
00127     for (reco::ConversionTrackCollection::const_iterator track=tC1.begin(); track!=tC1.end(); ++track){
00128       i++; 
00129       std::vector<const TrackingRecHit*>& iHits = rh1[track]; 
00130       unsigned nh1 = iHits.size();
00131       int j=-1;
00132       for (reco::ConversionTrackCollection::const_iterator track2=tC2.begin(); track2!=tC2.end(); ++track2){
00133         j++;
00134         std::vector<const TrackingRecHit*>& jHits = rh2[track2]; 
00135         unsigned nh2 = jHits.size();
00136         int noverlap=0;
00137         int firstoverlap=0;
00138         for ( unsigned ih=0; ih<nh1; ++ih ) { 
00139           const TrackingRecHit* it = iHits[ih];
00140           if (it->isValid()){
00141             int jj=-1;
00142             for ( unsigned jh=0; jh<nh2; ++jh ) { 
00143               const TrackingRecHit* jt = jHits[jh];
00144               jj++;
00145               if (jt->isValid()){           
00146                 if ( it->sharesInput(jt,TrackingRecHit::some) ) {
00147                   noverlap++;
00148                   if ( allowFirstHitShare && ( ih == 0 ) && ( jh == 0 ) ) firstoverlap=1;
00149                 }
00150               }
00151             }
00152           }
00153         }
00154         int nhit1 = track->track()->numberOfValidHits();
00155         int nhit2 = track2->track()->numberOfValidHits();
00156         //std::cout << " trk1 trk2 nhits1 nhits2 nover " << i << " " << j << " " << track->numberOfValidHits() << " "  << track2->numberOfValidHits() << " " << noverlap << " " << fi << " " << fj  <<std::endl;
00157         if ( (noverlap-firstoverlap) > (std::min(nhit1,nhit2)-firstoverlap)*shareFrac ) {
00158           if ( nhit1 > nhit2 ){
00159             selected2[j]=0; 
00160             //std::cout << " removing L2 trk in pair " << std::endl;
00161           }else{
00162             if ( nhit1 < nhit2 ){
00163               selected1[i]=0; 
00164               //std::cout << " removing L1 trk in pair " << std::endl;
00165             }else{
00166               //std::cout << " removing worst chisq in pair " << track->normalizedChi2() << " " << track2->normalizedChi2() << std::endl;
00167               if (track->track()->normalizedChi2() > track2->track()->normalizedChi2()) {
00168                 selected1[i]=0;
00169               }else {
00170                 selected2[j]=0;
00171               }
00172             }//end fi > or = fj
00173           }//end fi < fj
00174         }//end got a duplicate
00175       }//end track2 loop
00176     }//end track loop
00177    }//end more than 1 track
00178 
00179   //
00180   //  output selected tracks - if any
00181   //
00182    
00183    if ( 0<tC1.size() ){
00184      i=0;
00185      for (reco::ConversionTrackCollection::const_iterator track=tC1.begin(); track!=tC1.end(); 
00186           ++track, ++i){
00187       //don't store tracks rejected as duplicates
00188       if ( outputPreferCollection==0 || ( (outputPreferCollection==3 || outputPreferCollection==2) && !selected1[i]) ){
00189         continue;
00190       }
00191       //fill the TrackCollection
00192       outputTrks->push_back(*track);
00193       //clear flags for tracks rejected as duplicates
00194       if ( trackerOnlyPreferCollection==0 || ( (trackerOnlyPreferCollection==3 || trackerOnlyPreferCollection==2) && !selected1[i]) ){
00195         outputTrks->back().setIsTrackerOnly(false);
00196       }
00197       if ( arbitratedEcalSeededPreferCollection==0 || ( (arbitratedEcalSeededPreferCollection==3 || arbitratedEcalSeededPreferCollection==2) && !selected1[i]) ){
00198         outputTrks->back().setIsArbitratedEcalSeeded(false);
00199       }      
00200       if ( arbitratedMergedPreferCollection==0 || ( (arbitratedMergedPreferCollection==3 || arbitratedMergedPreferCollection==2) && !selected1[i]) ){
00201         outputTrks->back().setIsArbitratedMerged(false);
00202       }    
00203       if ( arbitratedMergedEcalGeneralPreferCollection==0 || ( (arbitratedMergedEcalGeneralPreferCollection==3 || arbitratedMergedEcalGeneralPreferCollection==2) && !selected1[i]) ){
00204         outputTrks->back().setIsArbitratedMergedEcalGeneral(false);
00205       }       
00206     }//end faux loop over tracks
00207    }//end more than 0 track
00208 
00209    //Fill the trajectories, etc. for 1st collection
00210  
00211    
00212    if ( 0<tC2.size() ){
00213     i=0;
00214     for (reco::ConversionTrackCollection::const_iterator track=tC2.begin(); track!=tC2.end();
00215          ++track, ++i){
00216       //don't store tracks rejected as duplicates
00217       if ( outputPreferCollection==0 || ( (outputPreferCollection==3 || outputPreferCollection==1) && !selected2[i]) ){
00218         continue;
00219       }
00220       //fill the TrackCollection
00221       outputTrks->push_back(*track);
00222       //clear flags for tracks rejected as duplicates
00223       if ( trackerOnlyPreferCollection==0 || ( (trackerOnlyPreferCollection==3 || trackerOnlyPreferCollection==1) && !selected2[i]) ){
00224         outputTrks->back().setIsTrackerOnly(false);
00225       }
00226       if ( arbitratedEcalSeededPreferCollection==0 || ( (arbitratedEcalSeededPreferCollection==3 || arbitratedEcalSeededPreferCollection==1) && !selected2[i]) ){
00227         outputTrks->back().setIsArbitratedEcalSeeded(false);
00228       }      
00229       if ( arbitratedMergedPreferCollection==0 || ( (arbitratedMergedPreferCollection==3 || arbitratedMergedPreferCollection==1) && !selected2[i]) ){
00230         outputTrks->back().setIsArbitratedMerged(false);
00231       }      
00232       if ( arbitratedMergedEcalGeneralPreferCollection==0 || ( (arbitratedMergedEcalGeneralPreferCollection==3 || arbitratedMergedEcalGeneralPreferCollection==1) && !selected2[i]) ){
00233         outputTrks->back().setIsArbitratedMergedEcalGeneral(false);
00234       }            
00235 
00236     }//end faux loop over tracks
00237    }//end more than 0 track
00238  
00239     e.put(outputTrks);
00240     return;
00241 
00242   }//end produce