CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoEgamma/EgammaElectronProducers/plugins/GsfElectronCoreBaseProducer.cc

Go to the documentation of this file.
00001 
00002 #include "GsfElectronCoreBaseProducer.h"
00003 
00004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00005 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
00006 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
00007 
00008 #include "DataFormats/EgammaCandidates/interface/GsfElectronCoreFwd.h"
00009 #include "DataFormats/EgammaCandidates/interface/GsfElectronCore.h"
00010 #include "DataFormats/ParticleFlowReco/interface/GsfPFRecTrack.h"
00011 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
00012 #include "DataFormats/EgammaReco/interface/ElectronSeedFwd.h"
00013 #include "DataFormats/EgammaReco/interface/ElectronSeed.h"
00014 //#include "DataFormats/Common/interface/ValueMap.h"
00015 
00016 //#include <map>
00017 
00018 using namespace reco ;
00019 
00020 void GsfElectronCoreBaseProducer::fillDescription( edm::ParameterSetDescription & desc )
00021  {
00022   desc.add<edm::InputTag>("gsfPfRecTracks",edm::InputTag("pfTrackElec")) ;
00023   desc.add<edm::InputTag>("gsfTracks",edm::InputTag("electronGsfTracks")) ;
00024   desc.add<edm::InputTag>("ctfTracks",edm::InputTag("generalTracks")) ;
00025   desc.add<bool>("useGsfPfRecTracks",true) ;
00026  }
00027 
00028 GsfElectronCoreBaseProducer::GsfElectronCoreBaseProducer( const edm::ParameterSet & config )
00029  {
00030   produces<GsfElectronCoreCollection>() ;
00031   gsfPfRecTracksTag_ = config.getParameter<edm::InputTag>("gsfPfRecTracks") ;
00032   gsfTracksTag_ = config.getParameter<edm::InputTag>("gsfTracks") ;
00033   ctfTracksTag_ = config.getParameter<edm::InputTag>("ctfTracks") ;
00034   useGsfPfRecTracks_ = config.getParameter<bool>("useGsfPfRecTracks") ;
00035  }
00036 
00037 GsfElectronCoreBaseProducer::~GsfElectronCoreBaseProducer()
00038  {}
00039 
00040 
00041 //=======================================================================================
00042 // For derived producers
00043 //=======================================================================================
00044 
00045 // to be called at the beginning of each new event
00046 void GsfElectronCoreBaseProducer::initEvent( edm::Event & event, const edm::EventSetup & setup )
00047  {
00048   if (useGsfPfRecTracks_)
00049    { event.getByLabel(gsfPfRecTracksTag_,gsfPfRecTracksH_) ; }
00050   event.getByLabel(gsfTracksTag_,gsfTracksH_) ;
00051   event.getByLabel(ctfTracksTag_,ctfTracksH_) ;
00052  }
00053 
00054 void GsfElectronCoreBaseProducer::fillElectronCore( reco::GsfElectronCore * eleCore )
00055  {
00056   const GsfTrackRef & gsfTrackRef = eleCore->gsfTrack() ;
00057 
00058   std::pair<TrackRef,float> ctfpair = getCtfTrackRef(gsfTrackRef) ;
00059   eleCore->setCtfTrack(ctfpair.first,ctfpair.second) ;
00060  }
00061 
00062 
00063 //=======================================================================================
00064 // Code from Puneeth Kalavase
00065 //=======================================================================================
00066 
00067 std::pair<TrackRef,float> GsfElectronCoreBaseProducer::getCtfTrackRef
00068  ( const GsfTrackRef & gsfTrackRef )
00069  {
00070   float maxFracShared = 0;
00071   TrackRef ctfTrackRef = TrackRef() ;
00072   const TrackCollection * ctfTrackCollection = ctfTracksH_.product() ;
00073 
00074   // get the Hit Pattern for the gsfTrack
00075   const HitPattern& gsfHitPattern = gsfTrackRef->hitPattern();
00076 
00077   unsigned int counter ;
00078   TrackCollection::const_iterator ctfTkIter ;
00079   for ( ctfTkIter = ctfTrackCollection->begin() , counter = 0 ;
00080         ctfTkIter != ctfTrackCollection->end() ; ctfTkIter++, counter++ )
00081    {
00082 
00083     double dEta = gsfTrackRef->eta() - ctfTkIter->eta();
00084     double dPhi = gsfTrackRef->phi() - ctfTkIter->phi();
00085     double pi = acos(-1.);
00086     if(std::abs(dPhi) > pi) dPhi = 2*pi - std::abs(dPhi);
00087 
00088     // dont want to look at every single track in the event!
00089     if(sqrt(dEta*dEta + dPhi*dPhi) > 0.3) continue;
00090 
00091     unsigned int shared = 0 ;
00092     int gsfHitCounter = 0 ;
00093     int numGsfInnerHits = 0 ;
00094     int numCtfInnerHits = 0 ;
00095     // get the CTF Track Hit Pattern
00096     const HitPattern& ctfHitPattern = ctfTkIter->hitPattern() ;
00097 
00098     trackingRecHit_iterator elHitsIt ;
00099     for ( elHitsIt = gsfTrackRef->recHitsBegin() ;
00100           elHitsIt != gsfTrackRef->recHitsEnd() ;
00101           elHitsIt++, gsfHitCounter++ )
00102      {
00103       if(!((**elHitsIt).isValid()))  //count only valid Hits
00104        { continue ; }
00105 
00106       // look only in the pixels/TIB/TID
00107       uint32_t gsfHit = gsfHitPattern.getHitPattern(gsfHitCounter) ;
00108       if (!(gsfHitPattern.pixelHitFilter(gsfHit) ||
00109           gsfHitPattern.stripTIBHitFilter(gsfHit) ||
00110           gsfHitPattern.stripTIDHitFilter(gsfHit) ) )
00111        { continue ; }
00112 
00113       numGsfInnerHits++ ;
00114 
00115       int ctfHitsCounter = 0 ;
00116       numCtfInnerHits = 0 ;
00117       trackingRecHit_iterator ctfHitsIt ;
00118       for ( ctfHitsIt = ctfTkIter->recHitsBegin() ;
00119             ctfHitsIt != ctfTkIter->recHitsEnd() ;
00120             ctfHitsIt++, ctfHitsCounter++ )
00121        {
00122         if(!((**ctfHitsIt).isValid())) //count only valid Hits!
00123          { continue ; }
00124 
00125       uint32_t ctfHit = ctfHitPattern.getHitPattern(ctfHitsCounter);
00126       if( !(ctfHitPattern.pixelHitFilter(ctfHit) ||
00127             ctfHitPattern.stripTIBHitFilter(ctfHit) ||
00128             ctfHitPattern.stripTIDHitFilter(ctfHit) ) )
00129        { continue ; }
00130 
00131       numCtfInnerHits++ ;
00132 
00133         if( (**elHitsIt).sharesInput(&(**ctfHitsIt),TrackingRecHit::all) )
00134          {
00135           shared++ ;
00136           break ;
00137          }
00138 
00139        } //ctfHits iterator
00140 
00141      } //gsfHits iterator
00142 
00143     if ((numGsfInnerHits==0)||(numCtfInnerHits==0))
00144      { continue ; }
00145 
00146     if ( static_cast<float>(shared)/std::min(numGsfInnerHits,numCtfInnerHits) > maxFracShared )
00147      {
00148       maxFracShared = static_cast<float>(shared)/std::min(numGsfInnerHits, numCtfInnerHits);
00149       ctfTrackRef = TrackRef(ctfTracksH_,counter);
00150      }
00151 
00152    } //ctfTrack iterator
00153 
00154   return make_pair(ctfTrackRef,maxFracShared) ;
00155  }
00156 
00157 
00158