File indexing completed on 2024-04-06 12:27:38
0001 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0002 #include "DQMServices/Core/interface/DQMStore.h"
0003 #include "FWCore/Framework/interface/Frameworkfwd.h"
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/Framework/interface/EventSetup.h"
0006 #include "FWCore/ServiceRegistry/interface/Service.h"
0007 #include "FWCore/Framework/interface/MakerMacros.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009 #include "DQMServices/Core/interface/DQMStore.h"
0010 #include "DataFormats/ParticleFlowReco/interface/PreIdFwd.h"
0011 #include "DataFormats/ParticleFlowReco/interface/PreId.h"
0012 #include "DataFormats/Common/interface/ValueMap.h"
0013
0014 class PreIdAnalyzer : public DQMEDAnalyzer {
0015 public:
0016 typedef dqm::legacy::DQMStore DQMStore;
0017 typedef dqm::legacy::MonitorElement MonitorElement;
0018
0019 explicit PreIdAnalyzer(const edm::ParameterSet&);
0020 ~PreIdAnalyzer() override = default;
0021
0022 void dqmBeginRun(const edm::Run&, const edm::EventSetup&) override {}
0023 void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
0024 void analyze(const edm::Event&, const edm::EventSetup&) override;
0025
0026 private:
0027 const edm::EDGetTokenT<edm::ValueMap<reco::PreIdRef> > preIdMapToken_;
0028 const edm::EDGetTokenT<reco::TrackCollection> trackToken_;
0029
0030 MonitorElement* tracksPt_;
0031 MonitorElement* tracksEta_;
0032 MonitorElement* tracksPtEcalMatch_;
0033 MonitorElement* tracksEtaEcalMatch_;
0034 MonitorElement* geomMatchChi2_;
0035 MonitorElement* geomMatchEop_;
0036
0037 MonitorElement* tracksChi2_;
0038 MonitorElement* tracksNhits_;
0039 MonitorElement* tracksPtFiltered_;
0040 MonitorElement* tracksEtaFiltered_;
0041 MonitorElement* tracksPtNotFiltered_;
0042 MonitorElement* tracksEtaNotFiltered_;
0043
0044 MonitorElement* tracksPtPreIded_;
0045 MonitorElement* tracksEtaPreIded_;
0046 MonitorElement* trackdpt_;
0047 MonitorElement* gsfChi2_;
0048 MonitorElement* chi2Ratio_;
0049 MonitorElement* mva_;
0050 };
0051
0052 PreIdAnalyzer::PreIdAnalyzer(const edm::ParameterSet& pset)
0053 : preIdMapToken_(consumes<edm::ValueMap<reco::PreIdRef> >(pset.getParameter<edm::InputTag>("PreIdMap"))),
0054 trackToken_(consumes<reco::TrackCollection>(pset.getParameter<edm::InputTag>("TrackCollection"))) {}
0055
0056 void PreIdAnalyzer::bookHistograms(DQMStore::IBooker& dbe, edm::Run const&, edm::EventSetup const&) {
0057
0058 tracksPt_ = dbe.book1D("TracksPt", "pT", 1000, 0, 100.);
0059 tracksEta_ = dbe.book1D("TracksEta", "eta", 50, -2.5, 2.5);
0060 tracksPtEcalMatch_ = dbe.book1D("TracksPtEcalMatch", "pT", 1000, 0, 100.);
0061 tracksEtaEcalMatch_ = dbe.book1D("TracksEtaEcalMatch", "eta", 50, -2.5, 2.5);
0062 tracksPtFiltered_ = dbe.book1D("TracksPtFiltered", "pT", 1000, 0, 100.);
0063 tracksEtaFiltered_ = dbe.book1D("TracksEtaFiltered", "eta", 50, -2.5, 2.5);
0064 tracksPtNotFiltered_ = dbe.book1D("TracksPtNotFiltered", "pT", 1000, 0, 100.);
0065 tracksEtaNotFiltered_ = dbe.book1D("TracksEtaNotFiltered", "eta", 50, -2.5, 2.5);
0066 tracksPtPreIded_ = dbe.book1D("TracksPtPreIded", "pT", 1000, 0, 100.);
0067 tracksEtaPreIded_ = dbe.book1D("TracksEtaPreIded", "eta", 50, -2.5, 2.5);
0068 tracksChi2_ = dbe.book1D("TracksChi2", "chi2", 100, 0, 10.);
0069 tracksNhits_ = dbe.book1D("TracksNhits", "Nhits", 30, -0.5, 29.5);
0070
0071 geomMatchChi2_ = dbe.book1D("geomMatchChi2", "Geom Chi2", 100, 0., 50.);
0072 geomMatchEop_ = dbe.book1D("geomMatchEop", "E/p", 100, 0., 5.);
0073 trackdpt_ = dbe.book1D("trackdpt", "dpt/pt", 100, 0., 5.);
0074 gsfChi2_ = dbe.book1D("gsfChi2", "GSF chi2", 100, 0., 10.);
0075 chi2Ratio_ = dbe.book1D("chi2Ratio", "Chi2 ratio", 100, 0., 10.);
0076 mva_ = dbe.book1D("mva", "mva", 100, -1., 1.);
0077 }
0078
0079 void PreIdAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0080 auto const& trackh = iEvent.getHandle(trackToken_);
0081 auto const& vmaph = iEvent.getHandle(preIdMapToken_);
0082
0083 const reco::TrackCollection& tracks = *(trackh.product());
0084 const edm::ValueMap<reco::PreIdRef>& preidMap = *(vmaph.product());
0085
0086 unsigned ntracks = tracks.size();
0087 for (unsigned itrack = 0; itrack < ntracks; ++itrack) {
0088 reco::TrackRef theTrackRef(trackh, itrack);
0089 tracksPt_->Fill(theTrackRef->pt());
0090 tracksEta_->Fill(theTrackRef->eta());
0091
0092 if (preidMap[theTrackRef].isNull())
0093 continue;
0094
0095 const reco::PreId& myPreId(*(preidMap[theTrackRef]));
0096 geomMatchChi2_->Fill(myPreId.geomMatching()[4]);
0097 geomMatchEop_->Fill(myPreId.eopMatch());
0098
0099 if (myPreId.ecalMatching()) {
0100 tracksPtEcalMatch_->Fill(theTrackRef->pt());
0101 tracksEtaEcalMatch_->Fill(theTrackRef->eta());
0102 } else {
0103 tracksChi2_->Fill(myPreId.kfChi2());
0104 tracksNhits_->Fill(myPreId.kfNHits());
0105 if (myPreId.trackFiltered()) {
0106 tracksPtFiltered_->Fill(theTrackRef->pt());
0107 tracksEtaFiltered_->Fill(theTrackRef->eta());
0108 trackdpt_->Fill(myPreId.dpt());
0109 gsfChi2_->Fill(myPreId.gsfChi2());
0110 chi2Ratio_->Fill(myPreId.chi2Ratio());
0111 mva_->Fill(myPreId.mva());
0112 } else {
0113 tracksPtNotFiltered_->Fill(theTrackRef->pt());
0114 tracksEtaNotFiltered_->Fill(theTrackRef->eta());
0115 }
0116 }
0117 if (myPreId.preIded()) {
0118 tracksPtPreIded_->Fill(theTrackRef->pt());
0119 tracksEtaPreIded_->Fill(theTrackRef->eta());
0120 }
0121 }
0122 }
0123
0124 DEFINE_FWK_MODULE(PreIdAnalyzer);