CMS 3D CMS Logo

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

Go to the documentation of this file.
00001 
00002 #include "GsfElectronCoreProducer.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/SuperClusterFwd.h"
00013 #include "DataFormats/EgammaReco/interface/ElectronSeedFwd.h"
00014 #include "DataFormats/EgammaReco/interface/ElectronSeed.h"
00015 #include "DataFormats/Common/interface/ValueMap.h"
00016 
00017 #include <map>
00018 
00019 using namespace reco ;
00020 
00021 // void GsfElectronCoreProducer::fillDescriptions( edm::ConfigurationDescriptions & descriptions )
00022 //  {
00023 //   edm::ParameterSetDescription desc ;
00024 //   GsfElectronCoreBaseProducer::fillDescriptions(desc) ;
00025 //   desc.add<edm::InputTag>("ecalDrivenGsfElectronCoresTag",edm::InputTag("ecalDrivenGsfElectronCores")) ;
00026 //   desc.add<edm::InputTag>("pflowDrivenGsfElectronCoresTag",edm::InputTag("pflowGsfElectronCores")) ;
00027 //   desc.add<edm::InputTag>("pfSuperClusters",edm::InputTag("pfElectronTranslator:pf")) ;
00028 //   desc.add<edm::InputTag>("pfSuperClusterTrackMap",edm::InputTag("pfElectronTranslator:pf")) ;
00029 //   descriptions.add("produceGsfElectronCores",desc) ;
00030 //  }
00031 
00032 GsfElectronCoreProducer::GsfElectronCoreProducer( const edm::ParameterSet & config )
00033  : GsfElectronCoreBaseProducer(config)
00034  {
00035   edCoresTag_ = config.getParameter<edm::InputTag>("ecalDrivenGsfElectronCoresTag") ;
00036   pfCoresTag_ = config.getParameter<edm::InputTag>("pflowGsfElectronCoresTag") ;
00037   pfSuperClustersTag_ = config.getParameter<edm::InputTag>("pfSuperClusters") ;
00038   pfSuperClusterTrackMapTag_ = config.getParameter<edm::InputTag>("pfSuperClusterTrackMap") ;
00039  }
00040 
00041 void GsfElectronCoreProducer::produce( edm::Event & event, const edm::EventSetup & setup )
00042  {
00043   // base input
00044   GsfElectronCoreBaseProducer::initEvent(event,setup) ;
00045 
00046   // transient output
00047   std::list<GsfElectronCore *> electrons ;
00048 
00049   // event input
00050   event.getByLabel(edCoresTag_,edCoresH_) ;
00051   event.getByLabel(pfCoresTag_,pfCoresH_) ;
00052   event.getByLabel(pfSuperClustersTag_,pfClustersH_) ;
00053   event.getByLabel(pfSuperClusterTrackMapTag_,pfClusterTracksH_) ;
00054 
00055   // loop on pure tracker driven tracks
00056   if (useGsfPfRecTracks_)
00057    {
00058     const GsfPFRecTrackCollection * gsfPfRecTrackCollection = gsfPfRecTracksH_.product() ;
00059     GsfPFRecTrackCollection::const_iterator gsfPfRecTrack ;
00060     for ( gsfPfRecTrack=gsfPfRecTrackCollection->begin() ;
00061           gsfPfRecTrack!=gsfPfRecTrackCollection->end() ;
00062           ++gsfPfRecTrack )
00063      {
00064       const GsfTrackRef gsfTrackRef = gsfPfRecTrack->gsfTrackRef() ;
00065       produceTrackerDrivenCore(gsfTrackRef,electrons) ;
00066      }
00067    }
00068   else
00069    {
00070     const GsfTrackCollection * gsfTrackCollection = gsfTracksH_.product() ;
00071     for ( unsigned int i=0 ; i<gsfTrackCollection->size() ; ++i )
00072      {
00073       const GsfTrackRef gsfTrackRef = edm::Ref<GsfTrackCollection>(gsfTracksH_,i) ;
00074       produceTrackerDrivenCore(gsfTrackRef,electrons) ;
00075      }
00076    }
00077 
00078   // clone ecal driven electrons
00079   const GsfElectronCoreCollection * edCoresCollection = edCoresH_.product() ;
00080   GsfElectronCoreCollection::const_iterator edCoreIter ;
00081   for
00082    ( edCoreIter = edCoresCollection->begin() ;
00083      edCoreIter != edCoresCollection->end() ;
00084      edCoreIter++ )
00085    { electrons.push_back(edCoreIter->clone()) ; }
00086 
00087   // add pflow info
00088   const GsfElectronCoreCollection * pfCoresCollection = pfCoresH_.product() ;
00089   GsfElectronCoreCollection::const_iterator pfCoreIter ;
00090   std::list<GsfElectronCore *>::iterator eleCore ;
00091   bool found ;
00092   for ( eleCore = electrons.begin() ; eleCore != electrons.end() ; eleCore++ )
00093    {
00094 //    (*eleCore)->setPflowSuperCluster((*pfClusterTracksH_)[(*eleCore)->gsfTrack()]) ;
00095     found = false ;
00096     for
00097      ( pfCoreIter = pfCoresCollection->begin() ;
00098        pfCoreIter != pfCoresCollection->end() ;
00099        pfCoreIter++ )
00100      {
00101       if (pfCoreIter->gsfTrack()==(*eleCore)->gsfTrack())
00102        {
00103         if (found)
00104          {
00105           edm::LogWarning("GsfElectronCoreProducer")<<"associated pfGsfElectronCore already found" ;
00106          }
00107         else
00108          {
00109           found = true ;
00110           (*eleCore)->setPflowSuperCluster(pfCoreIter->pflowSuperCluster()) ;
00111          }
00112        }
00113      }
00114    }
00115 
00116   // store
00117   std::auto_ptr<GsfElectronCoreCollection> collection(new GsfElectronCoreCollection) ;
00118   for ( eleCore = electrons.begin() ; eleCore != electrons.end() ; eleCore++ )
00119    {
00120     if ((*eleCore)->superCluster().isNull())
00121      { LogDebug("GsfElectronCoreProducer")<<"GsfTrack with no associated SuperCluster at all." ; }
00122     else
00123      { collection->push_back(**eleCore) ; }
00124     delete (*eleCore) ;
00125    }
00126   event.put(collection) ;
00127  }
00128 
00129 void GsfElectronCoreProducer::produceTrackerDrivenCore( const GsfTrackRef & gsfTrackRef, std::list<GsfElectronCore *> & electrons )
00130  {
00131   GsfElectronCore * eleCore = new GsfElectronCore(gsfTrackRef) ;
00132   if (eleCore->ecalDrivenSeed())
00133    { delete eleCore ; return ; }
00134   GsfElectronCoreBaseProducer::fillElectronCore(eleCore) ;
00135   electrons.push_back(eleCore) ;
00136  }
00137 
00138 GsfElectronCoreProducer::~GsfElectronCoreProducer()
00139  {}
00140