#include <ElectronGeneralAnalyzer.h>
Public Member Functions | |
virtual void | analyze (const edm::Event &e, const edm::EventSetup &c) |
virtual void | book () |
ElectronGeneralAnalyzer (const edm::ParameterSet &conf) | |
virtual | ~ElectronGeneralAnalyzer () |
Private Attributes | |
edm::InputTag | beamSpotTag_ |
edm::InputTag | electronCollection_ |
edm::InputTag | gsftrackCollection_ |
MonitorElement * | h1_ele_triggers |
MonitorElement * | h2_ele_beamSpotXvsY |
edm::InputTag | matchingObjectCollection_ |
MonitorElement * | py_ele_nClustersVsLs |
MonitorElement * | py_ele_nElectronsVsLs |
MonitorElement * | py_ele_nGsfTracksVsLs |
MonitorElement * | py_ele_nTracksVsLs |
MonitorElement * | py_ele_nVerticesVsLs |
edm::InputTag | trackCollection_ |
edm::InputTag | triggerResults_ |
edm::InputTag | vertexCollection_ |
Definition at line 19 of file ElectronGeneralAnalyzer.h.
ElectronGeneralAnalyzer::ElectronGeneralAnalyzer | ( | const edm::ParameterSet & | conf | ) | [explicit] |
Definition at line 29 of file ElectronGeneralAnalyzer.cc.
References beamSpotTag_, electronCollection_, edm::ParameterSet::getParameter(), gsftrackCollection_, matchingObjectCollection_, trackCollection_, triggerResults_, and vertexCollection_.
: ElectronDqmAnalyzerBase(conf) { // collection input tags electronCollection_ = conf.getParameter<edm::InputTag>("ElectronCollection"); matchingObjectCollection_ = conf.getParameter<edm::InputTag>("MatchingObjectCollection"); trackCollection_ = conf.getParameter<edm::InputTag>("TrackCollection"); vertexCollection_ = conf.getParameter<edm::InputTag>("VertexCollection"); gsftrackCollection_ = conf.getParameter<edm::InputTag>("GsfTrackCollection"); beamSpotTag_ = conf.getParameter<edm::InputTag>("BeamSpot"); triggerResults_ = conf.getParameter<edm::InputTag>("TriggerResults"); // // for trigger // HLTPathsByName_= conf.getParameter<std::vector<std::string > >("HltPaths"); // HLTPathsByIndex_.resize(HLTPathsByName_.size()); }
ElectronGeneralAnalyzer::~ElectronGeneralAnalyzer | ( | ) | [virtual] |
Definition at line 46 of file ElectronGeneralAnalyzer.cc.
{}
void ElectronGeneralAnalyzer::analyze | ( | const edm::Event & | e, |
const edm::EventSetup & | c | ||
) | [virtual] |
Reimplemented from ElectronDqmAnalyzerBase.
Definition at line 60 of file ElectronGeneralAnalyzer.cc.
References beamSpotTag_, electronCollection_, edm::EventID::event(), MonitorElement::Fill(), edm::Event::getByLabel(), gsfElectrons_cfi::gsfElectrons, gsftrackCollection_, h1_ele_triggers, h2_ele_beamSpotXvsY, i, edm::EventBase::id(), edm::HandleBase::isValid(), edm::EventBase::luminosityBlock(), matchingObjectCollection_, n, reco::BeamSpot::position(), edm::Handle< T >::product(), py_ele_nClustersVsLs, py_ele_nElectronsVsLs, py_ele_nGsfTracksVsLs, py_ele_nTracksVsLs, py_ele_nVerticesVsLs, edm::EventID::run(), trackCollection_, testEve_cfg::tracks, patRefSel_triggerSelection_cff::triggerResults, triggerResults_, and vertexCollection_.
{ edm::Handle<GsfElectronCollection> gsfElectrons ; iEvent.getByLabel(electronCollection_,gsfElectrons) ; edm::Handle<reco::SuperClusterCollection> recoClusters ; iEvent.getByLabel(matchingObjectCollection_,recoClusters) ; edm::Handle<reco::TrackCollection> tracks; iEvent.getByLabel(trackCollection_,tracks); edm::Handle<reco::GsfTrackCollection> gsfTracks; iEvent.getByLabel(gsftrackCollection_,gsfTracks); edm::Handle<reco::VertexCollection> vertices; iEvent.getByLabel(vertexCollection_,vertices); edm::Handle<reco::BeamSpot> recoBeamSpotHandle ; iEvent.getByLabel(beamSpotTag_,recoBeamSpotHandle) ; const BeamSpot bs = *recoBeamSpotHandle ; int ievt = iEvent.id().event(); int irun = iEvent.id().run(); int ils = iEvent.luminosityBlock(); edm::LogInfo("ElectronGeneralAnalyzer::analyze") <<"Treating "<<gsfElectrons.product()->size()<<" electrons" <<" from event "<<ievt<<" in run "<<irun<<" and lumiblock "<<ils ; h2_ele_beamSpotXvsY->Fill(bs.position().x(),bs.position().y()); py_ele_nElectronsVsLs->Fill(float(ils),(*gsfElectrons).size()); py_ele_nClustersVsLs->Fill(float(ils),(*recoClusters).size()); py_ele_nGsfTracksVsLs->Fill(float(ils),(*gsfTracks).size()); py_ele_nTracksVsLs->Fill(float(ils),(*tracks).size()); py_ele_nVerticesVsLs->Fill(float(ils),(*vertices).size()); // trigger edm::Handle<edm::TriggerResults> triggerResults ; iEvent.getByLabel(triggerResults_,triggerResults) ; if (triggerResults.isValid()) { unsigned int i, n = triggerResults->size() ; for ( i=0 ; i!=n ; ++i ) { if (triggerResults->accept(i)) { h1_ele_triggers->Fill(float(i)) ; } } } }
void ElectronGeneralAnalyzer::book | ( | ) | [virtual] |
Reimplemented from ElectronDqmAnalyzerBase.
Definition at line 49 of file ElectronGeneralAnalyzer.cc.
References ElectronDqmAnalyzerBase::bookH1(), ElectronDqmAnalyzerBase::bookH2(), ElectronDqmAnalyzerBase::bookP1(), h1_ele_triggers, h2_ele_beamSpotXvsY, py_ele_nClustersVsLs, py_ele_nElectronsVsLs, py_ele_nGsfTracksVsLs, py_ele_nTracksVsLs, and py_ele_nVerticesVsLs.
{ h2_ele_beamSpotXvsY = bookH2("beamSpotXvsY","beam spot x vs y",100,-0.2,0.2,100,-0.2,0.2,"x (cm)","y (cm)") ; py_ele_nElectronsVsLs = bookP1("nElectronsVsLs","# gsf electrons vs LS",150,0.,1500.,0.,20.,"LS","<N_{ele}>") ; py_ele_nClustersVsLs = bookP1("nClustersVsLs","# clusters vs LS",150,0.,1500.,0.,100.,"LS","<N_{SC}>") ; py_ele_nGsfTracksVsLs = bookP1("nGsfTracksVsLs","# gsf tracks vs LS",150,0.,1500.,0.,20.,"LS","<N_{GSF tk}>") ; py_ele_nTracksVsLs = bookP1("nTracksVsLs","# tracks vs LS",150,0.,1500.,0.,100.,"LS","<N_{gen tk}>") ; py_ele_nVerticesVsLs = bookP1("nVerticesVsLs","# vertices vs LS",150,0.,1500.,0.,10.,"LS","<N_{vert}>") ; h1_ele_triggers = bookH1("triggers","hlt triggers",256,0.,256.,"HLT bit") ; }
Definition at line 41 of file ElectronGeneralAnalyzer.h.
Referenced by analyze(), and ElectronGeneralAnalyzer().
Definition at line 36 of file ElectronGeneralAnalyzer.h.
Referenced by analyze(), and ElectronGeneralAnalyzer().
Definition at line 38 of file ElectronGeneralAnalyzer.h.
Referenced by analyze(), and ElectronGeneralAnalyzer().
Definition at line 64 of file ElectronGeneralAnalyzer.h.
Definition at line 58 of file ElectronGeneralAnalyzer.h.
Definition at line 37 of file ElectronGeneralAnalyzer.h.
Referenced by analyze(), and ElectronGeneralAnalyzer().
Definition at line 60 of file ElectronGeneralAnalyzer.h.
Definition at line 59 of file ElectronGeneralAnalyzer.h.
Definition at line 61 of file ElectronGeneralAnalyzer.h.
Definition at line 62 of file ElectronGeneralAnalyzer.h.
Definition at line 63 of file ElectronGeneralAnalyzer.h.
Definition at line 39 of file ElectronGeneralAnalyzer.h.
Referenced by analyze(), and ElectronGeneralAnalyzer().
Definition at line 44 of file ElectronGeneralAnalyzer.h.
Referenced by analyze(), and ElectronGeneralAnalyzer().
Definition at line 40 of file ElectronGeneralAnalyzer.h.
Referenced by analyze(), and ElectronGeneralAnalyzer().