File indexing completed on 2024-04-06 12:09:21
0001
0002 #include "DQMOffline/EGamma/plugins/ElectronGeneralAnalyzer.h"
0003
0004 #include "DQMServices/Core/interface/DQMStore.h"
0005
0006 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
0007 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0008 #include "DataFormats/EgammaReco/interface/BasicClusterFwd.h"
0009 #include "DataFormats/Common/interface/Handle.h"
0010 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0011 #include "DataFormats/Common/interface/TriggerResults.h"
0012
0013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0014 #include "FWCore/Framework/interface/Event.h"
0015 #include "FWCore/Framework/interface/MakerMacros.h"
0016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0017 #include "FWCore/ServiceRegistry/interface/Service.h"
0018
0019
0020
0021
0022 #include <iostream>
0023
0024 using namespace reco;
0025
0026 ElectronGeneralAnalyzer::ElectronGeneralAnalyzer(const edm::ParameterSet& conf) : ElectronDqmAnalyzerBase(conf) {
0027
0028 electronCollection_ = consumes<GsfElectronCollection>(conf.getParameter<edm::InputTag>("ElectronCollection"));
0029 matchingObjectCollection_ =
0030 consumes<reco::SuperClusterCollection>(conf.getParameter<edm::InputTag>("MatchingObjectCollection"));
0031 trackCollection_ = consumes<reco::TrackCollection>(conf.getParameter<edm::InputTag>("TrackCollection"));
0032 vertexCollection_ = consumes<reco::VertexCollection>(conf.getParameter<edm::InputTag>("VertexCollection"));
0033 gsftrackCollection_ = consumes<reco::GsfTrackCollection>(conf.getParameter<edm::InputTag>("GsfTrackCollection"));
0034 beamSpotTag_ = consumes<reco::BeamSpot>(conf.getParameter<edm::InputTag>("BeamSpot"));
0035 }
0036
0037 ElectronGeneralAnalyzer::~ElectronGeneralAnalyzer() {}
0038
0039 void ElectronGeneralAnalyzer::bookHistograms(DQMStore::IBooker& iBooker, edm::Run const&, edm::EventSetup const&) {
0040 h2_ele_beamSpotXvsY =
0041 bookH2(iBooker, "beamSpotXvsY", "beam spot x vs y", 100, -0.2, 0.2, 100, -0.2, 0.2, "x (cm)", "y (cm)");
0042 py_ele_nElectronsVsLs =
0043 bookP1(iBooker, "nElectronsVsLs", "# gsf electrons vs LS", 150, 0., 1500., 0., 20., "LS", "<N_{ele}>");
0044 py_ele_nClustersVsLs =
0045 bookP1(iBooker, "nClustersVsLs", "# clusters vs LS", 150, 0., 1500., 0., 100., "LS", "<N_{SC}>");
0046 py_ele_nGsfTracksVsLs =
0047 bookP1(iBooker, "nGsfTracksVsLs", "# gsf tracks vs LS", 150, 0., 1500., 0., 20., "LS", "<N_{GSF tk}>");
0048 py_ele_nTracksVsLs = bookP1(iBooker, "nTracksVsLs", "# tracks vs LS", 150, 0., 1500., 0., 100., "LS", "<N_{gen tk}>");
0049 py_ele_nVerticesVsLs =
0050 bookP1(iBooker, "nVerticesVsLs", "# vertices vs LS", 150, 0., 1500., 0., 10., "LS", "<N_{vert}>");
0051 }
0052
0053 void ElectronGeneralAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0054 edm::Handle<GsfElectronCollection> gsfElectrons;
0055 iEvent.getByToken(electronCollection_, gsfElectrons);
0056 edm::Handle<reco::SuperClusterCollection> recoClusters;
0057 iEvent.getByToken(matchingObjectCollection_, recoClusters);
0058 edm::Handle<reco::TrackCollection> tracks;
0059 iEvent.getByToken(trackCollection_, tracks);
0060 edm::Handle<reco::GsfTrackCollection> gsfTracks;
0061 iEvent.getByToken(gsftrackCollection_, gsfTracks);
0062 edm::Handle<reco::VertexCollection> vertices;
0063 iEvent.getByToken(vertexCollection_, vertices);
0064 edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
0065 iEvent.getByToken(beamSpotTag_, recoBeamSpotHandle);
0066 const BeamSpot bs = *recoBeamSpotHandle;
0067
0068 edm::EventNumber_t ievt = iEvent.id().event();
0069 edm::RunNumber_t irun = iEvent.id().run();
0070 edm::LuminosityBlockNumber_t ils = iEvent.luminosityBlock();
0071
0072 edm::LogInfo("ElectronGeneralAnalyzer::analyze")
0073 << "Treating " << gsfElectrons.product()->size() << " electrons"
0074 << " from event " << ievt << " in run " << irun << " and lumiblock " << ils;
0075
0076 h2_ele_beamSpotXvsY->Fill(bs.position().x(), bs.position().y());
0077 py_ele_nElectronsVsLs->Fill(float(ils), (*gsfElectrons).size());
0078 py_ele_nClustersVsLs->Fill(float(ils), (*recoClusters).size());
0079 py_ele_nGsfTracksVsLs->Fill(float(ils), (*gsfTracks).size());
0080 py_ele_nTracksVsLs->Fill(float(ils), (*tracks).size());
0081 py_ele_nVerticesVsLs->Fill(float(ils), (*vertices).size());
0082 }