#include <GsfElectronCoreProducer.h>
Public Member Functions | |
GsfElectronCoreProducer (const edm::ParameterSet &) | |
virtual void | produce (edm::Event &, const edm::EventSetup &) |
virtual | ~GsfElectronCoreProducer () |
Private Member Functions | |
void | produceTrackerDrivenCore (const reco::GsfTrackRef &gsfTrackRef, std::list< reco::GsfElectronCore * > &electrons) |
Private Attributes | |
edm::Handle < reco::GsfElectronCoreCollection > | edCoresH_ |
edm::InputTag | edCoresTag_ |
edm::Handle < reco::SuperClusterCollection > | pfClustersH_ |
edm::Handle< edm::ValueMap < reco::SuperClusterRef > > | pfClusterTracksH_ |
edm::Handle < reco::GsfElectronCoreCollection > | pfCoresH_ |
edm::InputTag | pfCoresTag_ |
edm::InputTag | pfSuperClustersTag_ |
edm::InputTag | pfSuperClusterTrackMapTag_ |
Definition at line 9 of file GsfElectronCoreProducer.h.
GsfElectronCoreProducer::GsfElectronCoreProducer | ( | const edm::ParameterSet & | config | ) | [explicit] |
Definition at line 32 of file GsfElectronCoreProducer.cc.
References edCoresTag_, edm::ParameterSet::getParameter(), pfCoresTag_, pfSuperClustersTag_, and pfSuperClusterTrackMapTag_.
: GsfElectronCoreBaseProducer(config) { edCoresTag_ = config.getParameter<edm::InputTag>("ecalDrivenGsfElectronCoresTag") ; pfCoresTag_ = config.getParameter<edm::InputTag>("pflowGsfElectronCoresTag") ; pfSuperClustersTag_ = config.getParameter<edm::InputTag>("pfSuperClusters") ; pfSuperClusterTrackMapTag_ = config.getParameter<edm::InputTag>("pfSuperClusterTrackMap") ; }
GsfElectronCoreProducer::~GsfElectronCoreProducer | ( | ) | [virtual] |
Definition at line 138 of file GsfElectronCoreProducer.cc.
{}
void GsfElectronCoreProducer::produce | ( | edm::Event & | event, |
const edm::EventSetup & | setup | ||
) | [virtual] |
Implements edm::EDProducer.
Definition at line 41 of file GsfElectronCoreProducer.cc.
References runEdmFileComparison::collection, edCoresH_, edCoresTag_, HI_PhotonSkim_cff::electrons, funct::false, newFWLiteAna::found, GsfElectronCoreBaseProducer::gsfPfRecTracksH_, GsfElectronCoreBaseProducer::gsfTracksH_, i, GsfElectronCoreBaseProducer::initEvent(), LogDebug, pfClustersH_, pfClusterTracksH_, pfCoresH_, pfCoresTag_, pfSuperClustersTag_, pfSuperClusterTrackMapTag_, produceTrackerDrivenCore(), edm::Handle< T >::product(), funct::true, and GsfElectronCoreBaseProducer::useGsfPfRecTracks_.
{ // base input GsfElectronCoreBaseProducer::initEvent(event,setup) ; // transient output std::list<GsfElectronCore *> electrons ; // event input event.getByLabel(edCoresTag_,edCoresH_) ; event.getByLabel(pfCoresTag_,pfCoresH_) ; event.getByLabel(pfSuperClustersTag_,pfClustersH_) ; event.getByLabel(pfSuperClusterTrackMapTag_,pfClusterTracksH_) ; // loop on pure tracker driven tracks if (useGsfPfRecTracks_) { const GsfPFRecTrackCollection * gsfPfRecTrackCollection = gsfPfRecTracksH_.product() ; GsfPFRecTrackCollection::const_iterator gsfPfRecTrack ; for ( gsfPfRecTrack=gsfPfRecTrackCollection->begin() ; gsfPfRecTrack!=gsfPfRecTrackCollection->end() ; ++gsfPfRecTrack ) { const GsfTrackRef gsfTrackRef = gsfPfRecTrack->gsfTrackRef() ; produceTrackerDrivenCore(gsfTrackRef,electrons) ; } } else { const GsfTrackCollection * gsfTrackCollection = gsfTracksH_.product() ; for ( unsigned int i=0 ; i<gsfTrackCollection->size() ; ++i ) { const GsfTrackRef gsfTrackRef = edm::Ref<GsfTrackCollection>(gsfTracksH_,i) ; produceTrackerDrivenCore(gsfTrackRef,electrons) ; } } // clone ecal driven electrons const GsfElectronCoreCollection * edCoresCollection = edCoresH_.product() ; GsfElectronCoreCollection::const_iterator edCoreIter ; for ( edCoreIter = edCoresCollection->begin() ; edCoreIter != edCoresCollection->end() ; edCoreIter++ ) { electrons.push_back(edCoreIter->clone()) ; } // add pflow info const GsfElectronCoreCollection * pfCoresCollection = pfCoresH_.product() ; GsfElectronCoreCollection::const_iterator pfCoreIter ; std::list<GsfElectronCore *>::iterator eleCore ; bool found ; for ( eleCore = electrons.begin() ; eleCore != electrons.end() ; eleCore++ ) { // (*eleCore)->setPflowSuperCluster((*pfClusterTracksH_)[(*eleCore)->gsfTrack()]) ; found = false ; for ( pfCoreIter = pfCoresCollection->begin() ; pfCoreIter != pfCoresCollection->end() ; pfCoreIter++ ) { if (pfCoreIter->gsfTrack()==(*eleCore)->gsfTrack()) { if (found) { edm::LogWarning("GsfElectronCoreProducer")<<"associated pfGsfElectronCore already found" ; } else { found = true ; (*eleCore)->setPflowSuperCluster(pfCoreIter->pflowSuperCluster()) ; } } } } // store std::auto_ptr<GsfElectronCoreCollection> collection(new GsfElectronCoreCollection) ; for ( eleCore = electrons.begin() ; eleCore != electrons.end() ; eleCore++ ) { if ((*eleCore)->superCluster().isNull()) { LogDebug("GsfElectronCoreProducer")<<"GsfTrack with no associated SuperCluster at all." ; } else { collection->push_back(**eleCore) ; } delete (*eleCore) ; } event.put(collection) ; }
void GsfElectronCoreProducer::produceTrackerDrivenCore | ( | const reco::GsfTrackRef & | gsfTrackRef, |
std::list< reco::GsfElectronCore * > & | electrons | ||
) | [private] |
Definition at line 129 of file GsfElectronCoreProducer.cc.
References reco::GsfElectronCore::ecalDrivenSeed(), and GsfElectronCoreBaseProducer::fillElectronCore().
Referenced by produce().
{ GsfElectronCore * eleCore = new GsfElectronCore(gsfTrackRef) ; if (eleCore->ecalDrivenSeed()) { delete eleCore ; return ; } GsfElectronCoreBaseProducer::fillElectronCore(eleCore) ; electrons.push_back(eleCore) ; }
Definition at line 26 of file GsfElectronCoreProducer.h.
Referenced by produce().
Definition at line 21 of file GsfElectronCoreProducer.h.
Referenced by GsfElectronCoreProducer(), and produce().
Definition at line 28 of file GsfElectronCoreProducer.h.
Referenced by produce().
edm::Handle<edm::ValueMap<reco::SuperClusterRef> > GsfElectronCoreProducer::pfClusterTracksH_ [private] |
Definition at line 29 of file GsfElectronCoreProducer.h.
Referenced by produce().
Definition at line 27 of file GsfElectronCoreProducer.h.
Referenced by produce().
Definition at line 22 of file GsfElectronCoreProducer.h.
Referenced by GsfElectronCoreProducer(), and produce().
Definition at line 23 of file GsfElectronCoreProducer.h.
Referenced by GsfElectronCoreProducer(), and produce().
Definition at line 24 of file GsfElectronCoreProducer.h.
Referenced by GsfElectronCoreProducer(), and produce().