Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 10:58:43

0001 /*
0002   TagAndProbeBtagTriggerMonitor DQM code
0003 */
0004 //
0005 // Originally created by:  Roberval Walsh
0006 //                         June 2017
0007 #include <memory>
0008 
0009 #include "CommonTools/TriggerUtils/interface/GenericTriggerEventFlag.h"
0010 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0011 #include "DQMServices/Core/interface/DQMStore.h"
0012 #include "DataFormats/BTauReco/interface/JetTag.h"
0013 #include "DataFormats/Common/interface/TriggerResults.h"
0014 #include "DataFormats/HLTReco/interface/TriggerEvent.h"
0015 #include "DataFormats/HLTReco/interface/TriggerObject.h"
0016 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
0017 #include "DataFormats/Math/interface/deltaR.h"
0018 #include "FWCore/Framework/interface/Frameworkfwd.h"
0019 #include "FWCore/Framework/interface/MakerMacros.h"
0020 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0021 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0022 #include "FWCore/Utilities/interface/transform.h"
0023 
0024 class TagAndProbeBtagTriggerMonitor : public DQMEDAnalyzer {
0025 public:
0026   TagAndProbeBtagTriggerMonitor(const edm::ParameterSet&);
0027   ~TagAndProbeBtagTriggerMonitor() throw() override;
0028 
0029 protected:
0030   void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
0031   void analyze(edm::Event const& iEvent, edm::EventSetup const& iSetup) override;
0032 
0033   struct JetRefCompare {
0034     inline bool operator()(const edm::RefToBase<reco::Jet>& j1, const edm::RefToBase<reco::Jet>& j2) const {
0035       return (j1.id() < j2.id()) || ((j1.id() == j2.id()) && (j1.key() < j2.key()));
0036     }
0037   };
0038   typedef std::map<edm::RefToBase<reco::Jet>, float, JetRefCompare> JetTagMap;
0039 
0040 private:
0041   const std::string folderName_;
0042 
0043   const bool requireValidHLTPaths_;
0044   bool hltPathsAreValid_;
0045 
0046   std::string processname_;
0047   std::string triggerobjbtag_;
0048 
0049   double jetPtmin_;
0050   double jetEtamax_;
0051   double tagBtagmin_;
0052   double probeBtagmin_;
0053 
0054   std::vector<double> jetPtbins_;
0055   std::vector<double> jetEtabins_;
0056   std::vector<double> jetPhibins_;
0057   std::vector<double> jetBtagbins_;
0058 
0059   edm::InputTag triggerSummaryLabel_;
0060 
0061   std::vector<edm::EDGetTokenT<reco::JetTagCollection> > jetTagTokens_;
0062 
0063   edm::EDGetTokenT<trigger::TriggerEvent> triggerSummaryToken_;
0064 
0065   MonitorElement* pt_jet1_;
0066   MonitorElement* pt_jet2_;
0067   MonitorElement* eta_jet1_;
0068   MonitorElement* eta_jet2_;
0069   MonitorElement* phi_jet1_;
0070   MonitorElement* phi_jet2_;
0071   MonitorElement* eta_phi_jet1_;
0072   MonitorElement* eta_phi_jet2_;
0073 
0074   MonitorElement* pt_probe_;
0075   MonitorElement* pt_probe_match_;
0076   MonitorElement* eta_probe_;
0077   MonitorElement* eta_probe_match_;
0078   MonitorElement* phi_probe_;
0079   MonitorElement* phi_probe_match_;
0080   MonitorElement* eta_phi_probe_;
0081   MonitorElement* eta_phi_probe_match_;
0082 
0083   MonitorElement* discr_offline_btag_jet1_;
0084   MonitorElement* discr_offline_btag_jet2_;
0085 
0086   std::unique_ptr<GenericTriggerEventFlag> genTriggerEventFlag_;  // tag & probe: trigger flag for num and den
0087 };
0088 
0089 TagAndProbeBtagTriggerMonitor::TagAndProbeBtagTriggerMonitor(const edm::ParameterSet& iConfig)
0090     : folderName_(iConfig.getParameter<std::string>("dirname")),
0091       requireValidHLTPaths_(iConfig.getParameter<bool>("requireValidHLTPaths")),
0092       hltPathsAreValid_(false) {
0093   processname_ = iConfig.getParameter<std::string>("processname");
0094   triggerobjbtag_ = iConfig.getParameter<std::string>("triggerobjbtag");
0095   jetPtmin_ = iConfig.getParameter<double>("jetPtMin");
0096   jetEtamax_ = iConfig.getParameter<double>("jetEtaMax");
0097   tagBtagmin_ = iConfig.getParameter<double>("tagBtagMin");
0098   probeBtagmin_ = iConfig.getParameter<double>("probeBtagMin");
0099   triggerSummaryLabel_ = iConfig.getParameter<edm::InputTag>("triggerSummary");
0100   triggerSummaryToken_ = consumes<trigger::TriggerEvent>(triggerSummaryLabel_);
0101   // New
0102   jetTagTokens_ =
0103       edm::vector_transform(iConfig.getParameter<std::vector<edm::InputTag> >("btagAlgos"),
0104                             [this](edm::InputTag const& tag) { return mayConsume<reco::JetTagCollection>(tag); });
0105 
0106   genTriggerEventFlag_ = std::make_unique<GenericTriggerEventFlag>(
0107       iConfig.getParameter<edm::ParameterSet>("genericTriggerEventPSet"), consumesCollector(), *this);
0108 
0109   jetPtbins_ = iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("jetPt");
0110   jetEtabins_ = iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("jetEta");
0111   jetPhibins_ = iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("jetPhi");
0112   jetBtagbins_ = iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("jetBtag");
0113 }
0114 
0115 TagAndProbeBtagTriggerMonitor::~TagAndProbeBtagTriggerMonitor() throw() {
0116   if (genTriggerEventFlag_) {
0117     genTriggerEventFlag_.reset();
0118   }
0119 }
0120 
0121 void TagAndProbeBtagTriggerMonitor::bookHistograms(DQMStore::IBooker& ibooker,
0122                                                    edm::Run const& iRun,
0123                                                    edm::EventSetup const& iSetup) {
0124   // Initialize the GenericTriggerEventFlag
0125   if (genTriggerEventFlag_ && genTriggerEventFlag_->on())
0126     genTriggerEventFlag_->initRun(iRun, iSetup);
0127 
0128   // check if every HLT path specified in numerator and denominator has a valid match in the HLT Menu
0129   hltPathsAreValid_ =
0130       (genTriggerEventFlag_ && genTriggerEventFlag_->on() && genTriggerEventFlag_->allHLTPathsAreValid());
0131 
0132   // if valid HLT paths are required,
0133   // create DQM outputs only if all paths are valid
0134   if (requireValidHLTPaths_ and (not hltPathsAreValid_)) {
0135     return;
0136   }
0137 
0138   std::string currentFolder = folderName_;
0139   ibooker.setCurrentFolder(currentFolder);
0140 
0141   int ptnbins = jetPtbins_.size() - 1;
0142   int etanbins = jetEtabins_.size() - 1;
0143   int phinbins = jetPhibins_.size() - 1;
0144   int btagnbins = jetBtagbins_.size() - 1;
0145 
0146   std::vector<float> fptbins(jetPtbins_.begin(), jetPtbins_.end());
0147   std::vector<float> fetabins(jetEtabins_.begin(), jetEtabins_.end());
0148   std::vector<float> fphibins(jetPhibins_.begin(), jetPhibins_.end());
0149   std::vector<float> fbtagbins(jetBtagbins_.begin(), jetBtagbins_.end());
0150 
0151   float* ptbins = &fptbins[0];
0152   float* etabins = &fetabins[0];
0153   float* phibins = &fphibins[0];
0154   float* btagbins = &fbtagbins[0];
0155 
0156   pt_jet1_ = ibooker.book1D("pt_jet1", "pt_jet1", ptnbins, ptbins);
0157   pt_jet2_ = ibooker.book1D("pt_jet2", "pt_jet2", ptnbins, ptbins);
0158   eta_jet1_ = ibooker.book1D("eta_jet1", "eta_jet1", etanbins, etabins);
0159   eta_jet2_ = ibooker.book1D("eta_jet2", "eta_jet2", etanbins, etabins);
0160   phi_jet1_ = ibooker.book1D("phi_jet1", "phi_jet1", phinbins, phibins);
0161   phi_jet2_ = ibooker.book1D("phi_jet2", "phi_jet2", phinbins, phibins);
0162   eta_phi_jet1_ = ibooker.book2D("eta_phi_jet1", "eta_phi_jet1", etanbins, etabins, phinbins, phibins);
0163   eta_phi_jet2_ = ibooker.book2D("eta_phi_jet2", "eta_phi_jet2", etanbins, etabins, phinbins, phibins);
0164 
0165   pt_probe_ = ibooker.book1D("pt_probe", "pt_probe", ptnbins, ptbins);
0166   pt_probe_match_ = ibooker.book1D("pt_probe_match", "pt_probe_match", ptnbins, ptbins);
0167   eta_probe_ = ibooker.book1D("eta_probe", "eta_probe", etanbins, etabins);
0168   eta_probe_match_ = ibooker.book1D("eta_probe_match", "eta_probe_match", etanbins, etabins);
0169   phi_probe_ = ibooker.book1D("phi_probe", "phi_probe", phinbins, phibins);
0170   phi_probe_match_ = ibooker.book1D("phi_probe_match", "phi_probe_match", phinbins, phibins);
0171   eta_phi_probe_ = ibooker.book2D("eta_phi_probe", "eta_phi_probe", etanbins, etabins, phinbins, phibins);
0172   eta_phi_probe_match_ = ibooker.book2D("eta_phi_probe_match", "eta_phi_match", etanbins, etabins, phinbins, phibins);
0173 
0174   discr_offline_btag_jet1_ = ibooker.book1D("discr_offline_btag_jet1", "discr_offline_btag_jet1", btagnbins, btagbins);
0175   discr_offline_btag_jet2_ = ibooker.book1D("discr_offline_btag_jet2", "discr_offline_btag_jet2", btagnbins, btagbins);
0176 }
0177 
0178 void TagAndProbeBtagTriggerMonitor::analyze(edm::Event const& iEvent, edm::EventSetup const& iSetup) {
0179   // if valid HLT paths are required,
0180   // analyze event only if all paths are valid
0181   if (requireValidHLTPaths_ and (not hltPathsAreValid_)) {
0182     return;
0183   }
0184 
0185   bool match1 = false;
0186   bool match2 = false;
0187 
0188   JetTagMap allJetBTagVals;
0189   for (const auto& jetTagToken : jetTagTokens_) {
0190     edm::Handle<reco::JetTagCollection> bjetHandle;
0191     iEvent.getByToken(jetTagToken, bjetHandle);
0192     if (not bjetHandle.isValid()) {
0193       edm::LogWarning("TagAndProbeBtagTriggerMonitor") << "B-Jet handle not valid, will skip event \n";
0194       return;
0195     }
0196 
0197     const reco::JetTagCollection& bTags = *(bjetHandle.product());
0198     for (const auto& i_jetTag : bTags) {
0199       const auto& jetRef = i_jetTag.first;
0200       const auto btagVal = i_jetTag.second;
0201       if (not std::isfinite(btagVal)) {
0202         continue;
0203       }
0204       if (allJetBTagVals.find(jetRef) != allJetBTagVals.end()) {
0205         allJetBTagVals.at(jetRef) += btagVal;
0206       } else {
0207         allJetBTagVals.insert(JetTagMap::value_type(jetRef, btagVal));
0208       }
0209     }
0210   }
0211 
0212   // applying selection for event; tag & probe -> selection  for all events
0213   if (genTriggerEventFlag_->on() && genTriggerEventFlag_->accept(iEvent, iSetup)) {
0214     if (allJetBTagVals.size() > 1) {
0215       auto jetBTagVal = allJetBTagVals.begin();
0216       auto jet1 = *dynamic_cast<const reco::Jet*>(jetBTagVal->first.get());
0217       auto btag1 = jetBTagVal->second;
0218 
0219       ++jetBTagVal;
0220       auto jet2 = *dynamic_cast<const reco::Jet*>(jetBTagVal->first.get());
0221       auto btag2 = jetBTagVal->second;
0222 
0223       if (jet1.pt() > jetPtmin_ && jet2.pt() > jetPtmin_ && fabs(jet1.eta()) < jetEtamax_ &&
0224           fabs(jet2.eta()) < jetEtamax_) {
0225         if (btag1 > tagBtagmin_ && btag2 > probeBtagmin_) {
0226           pt_jet1_->Fill(jet1.pt());
0227           pt_jet2_->Fill(jet2.pt());
0228           eta_jet1_->Fill(jet1.eta());
0229           eta_jet2_->Fill(jet2.eta());
0230           phi_jet1_->Fill(jet1.phi());
0231           phi_jet2_->Fill(jet2.phi());
0232           eta_phi_jet1_->Fill(jet1.eta(), jet1.phi());
0233           eta_phi_jet2_->Fill(jet2.eta(), jet2.phi());
0234           discr_offline_btag_jet1_->Fill(btag1);
0235           discr_offline_btag_jet2_->Fill(btag2);
0236 
0237           // trigger objects matching
0238           std::vector<trigger::TriggerObject> onlinebtags;
0239           edm::Handle<trigger::TriggerEvent> triggerEventHandler;
0240           iEvent.getByToken(triggerSummaryToken_, triggerEventHandler);
0241           const unsigned int filterIndex(
0242               triggerEventHandler->filterIndex(edm::InputTag(triggerobjbtag_, "", processname_)));
0243           if (filterIndex < triggerEventHandler->sizeFilters()) {
0244             const trigger::Keys& keys(triggerEventHandler->filterKeys(filterIndex));
0245             const trigger::TriggerObjectCollection& triggerObjects = triggerEventHandler->getObjects();
0246             for (auto& key : keys) {
0247               onlinebtags.reserve(onlinebtags.size() + keys.size());
0248               onlinebtags.push_back(triggerObjects[key]);
0249             }
0250           }
0251           for (auto const& to : onlinebtags) {
0252             if (reco::deltaR2(jet1, to) < 0.09)
0253               match1 = true;
0254             if (reco::deltaR2(jet2, to) < 0.09)
0255               match2 = true;
0256           }
0257 
0258           if (match1)  // jet1 is the tag
0259           {
0260             pt_probe_->Fill(jet2.pt());
0261             eta_probe_->Fill(jet2.eta());
0262             phi_probe_->Fill(jet2.phi());
0263             eta_phi_probe_->Fill(jet2.eta(), jet2.phi());
0264             if (match2)  // jet2 is the probe
0265             {
0266               pt_probe_match_->Fill(jet2.pt());
0267               eta_probe_match_->Fill(jet2.eta());
0268               phi_probe_match_->Fill(jet2.phi());
0269               eta_phi_probe_match_->Fill(jet2.eta(), jet2.phi());
0270             }
0271           }
0272         }
0273       }
0274     }
0275   }
0276 }
0277 
0278 DEFINE_FWK_MODULE(TagAndProbeBtagTriggerMonitor);