Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:07:37

0001 // system include files
0002 #include <memory>
0003 #include <sys/time.h>
0004 #include <cstdlib>
0005 #include <unordered_map>
0006 
0007 // user include files
0008 #include "DataFormats/Common/interface/TriggerResults.h"
0009 #include "DataFormats/HLTReco/interface/TriggerObject.h"
0010 #include "DataFormats/HLTReco/interface/TriggerEvent.h"
0011 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
0012 
0013 #include "FWCore/Utilities/interface/EDGetToken.h"
0014 #include "FWCore/Framework/interface/Frameworkfwd.h"
0015 #include "FWCore/Framework/interface/Event.h"
0016 #include "FWCore/Framework/interface/MakerMacros.h"
0017 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0018 #include "FWCore/ServiceRegistry/interface/Service.h"
0019 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0020 #include <DQMServices/Core/interface/DQMEDAnalyzer.h>
0021 
0022 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
0023 
0024 //for collections
0025 #include "HLTrigger/JetMET/interface/AlphaT.h"
0026 
0027 #include "DataFormats/Math/interface/deltaR.h"
0028 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0029 #include "DataFormats/BTauReco/interface/JetTag.h"
0030 #include "DataFormats/METReco/interface/MET.h"
0031 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
0032 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidateFwd.h"
0033 #include "DataFormats/JetReco/interface/PFJet.h"
0034 #include "DataFormats/JetReco/interface/CaloJet.h"
0035 
0036 #include "DQMServices/Core/interface/DQMStore.h"
0037 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0038 
0039 #include "TLorentzVector.h"
0040 
0041 struct hltPlot {
0042   typedef dqm::reco::MonitorElement MonitorElement;
0043 
0044   std::pair<MonitorElement*, bool> nME;
0045   std::pair<MonitorElement*, bool> etaME;
0046   std::pair<MonitorElement*, bool> phiME;
0047   std::pair<MonitorElement*, bool> ptME;
0048   std::pair<MonitorElement*, bool> massME;
0049   std::pair<MonitorElement*, bool> energyME;
0050   std::pair<MonitorElement*, bool> csvME;
0051   std::pair<MonitorElement*, bool> etaVSphiME;
0052   std::pair<MonitorElement*, bool> ptMEhep17;
0053   std::pair<MonitorElement*, bool> ptMEhem17;  // in harvesting step ratio
0054   std::pair<MonitorElement*, bool> mrME;
0055   std::pair<MonitorElement*, bool> rsqME;
0056   std::pair<MonitorElement*, bool> dxyME;
0057   std::pair<MonitorElement*, bool> dzME;
0058   std::pair<MonitorElement*, bool> dimassME;
0059   std::pair<MonitorElement*, bool> dRME;
0060   std::pair<MonitorElement*, bool> dRetaVSphiME;
0061   std::pair<MonitorElement*, bool> q1q2ME;
0062 
0063   std::string label;
0064   std::string pathNAME;
0065   int pathIDX;
0066   std::string moduleNAME;
0067 
0068   std::string xTITLE;
0069   std::vector<double> etaBINNING;
0070   std::vector<double> ptBINNING;
0071   std::vector<double> phiBINNING;
0072   std::vector<double> massBINNING;
0073   std::vector<double> dxyBINNING;
0074   std::vector<double> dzBINNING;
0075   std::vector<double> dimassBINNING;
0076 
0077   bool doPlot2D;
0078   bool doPlotETA;
0079   bool doPlotMASS;
0080   bool doPlotENERGY;
0081   bool doPlotHEP17;
0082   bool doPlotCSV;
0083   bool doCALO;
0084   bool doPF;
0085   bool doPlotRazor;
0086   bool doPlotDXY;
0087   bool doPlotDZ;
0088   bool doPlotDiMass;
0089 };
0090 //
0091 // class declaration
0092 //
0093 
0094 //using namespace edm;
0095 using std::unordered_map;
0096 
0097 class HLTObjectsMonitor : public DQMEDAnalyzer {
0098 public:
0099   explicit HLTObjectsMonitor(const edm::ParameterSet&);
0100   ~HLTObjectsMonitor() override = default;
0101 
0102   //      static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0103 
0104 private:
0105   void analyze(const edm::Event&, const edm::EventSetup&) override;
0106   void bookHistograms(DQMStore::IBooker& i, edm::Run const&, edm::EventSetup const&) override;
0107   void dqmBeginRun(edm::Run const&, edm::EventSetup const&) override;
0108 
0109   static hltPlot getPlotPSet(edm::ParameterSet pset);
0110   void getPSet();
0111   bool isHEP17(double eta, double phi);
0112   bool isHEM17(double eta, double phi);
0113 
0114   double dxyFinder(
0115       double, double, edm::Handle<reco::RecoChargedCandidateCollection>, edm::Handle<reco::BeamSpot>, double);
0116   double dzFinder(double, double, double, double, edm::Handle<reco::RecoChargedCandidateCollection>, double);
0117   // ----------member data ---------------------------
0118 
0119   std::string TopFolder_;
0120   std::string label_;
0121   std::string processName_;
0122   std::vector<edm::ParameterSet> plotPSETS_;
0123 
0124   HLTConfigProvider hltConfig_;
0125 
0126   MonitorElement* eventsPlot_;
0127   std::vector<hltPlot> hltPlots_;
0128   std::unordered_map<std::string, int> plotNamesToBins_;
0129 
0130   bool debug_;
0131 
0132   std::string mainFolder_;
0133   std::string backupFolder_;
0134 
0135   edm::EDGetTokenT<edm::TriggerResults> triggerResultsToken_;
0136   edm::EDGetTokenT<trigger::TriggerEvent> triggerEventToken_;
0137 
0138   edm::InputTag beamSpot_;
0139   edm::EDGetTokenT<reco::BeamSpot> beamSpotToken_;
0140 
0141   edm::EDGetTokenT<reco::JetTagCollection> caloJetBTagsToken_;
0142   edm::EDGetTokenT<reco::JetTagCollection> pfJetBTagsToken_;
0143 
0144   edm::InputTag muCandidates_;
0145   edm::EDGetTokenT<std::vector<reco::RecoChargedCandidate>> muCandidatesToken_;
0146   edm::InputTag eleCandidates_;
0147   edm::EDGetTokenT<std::vector<reco::RecoChargedCandidate>> eleCandidatesToken_;
0148 
0149   const double MASS_MU = .105658;
0150 
0151   struct MEbinning {
0152     int nbins;
0153     double xmin;
0154     double xmax;
0155   };
0156 
0157   double MAX_PHI = 3.2;
0158   int N_PHI = 64;
0159   const MEbinning phi_binning_{N_PHI, -MAX_PHI, MAX_PHI};
0160 
0161   double MAX_CSV = 1.;
0162   int N_CSV = 20;
0163   const MEbinning csv_binning_{N_CSV, -MAX_CSV, MAX_CSV};
0164 
0165   std::vector<double> phi_variable_binning_;
0166 
0167   /*
0168   HEP17 covers 
0169   - phi between 310° and 330° (-50° to -30°, or -0.87 t.52 rad)
0170   - eta between +1.3 and +3.0 (positive side only)
0171 */
0172   double MAX_PHI_HEP17 = -0.52;
0173   double MIN_PHI_HEP17 = -0.87;
0174   int N_PHI_HEP17 = 7;
0175   const MEbinning phi_binning_hep17_{N_PHI_HEP17, MIN_PHI_HEP17, MAX_PHI_HEP17};
0176   double MAX_ETA_HEP17 = 3.0;
0177   double MIN_ETA_HEP17 = 1.3;
0178   int N_ETA_HEP17 = 6;
0179   const MEbinning eta_binning_hep17_{N_ETA_HEP17, MIN_ETA_HEP17, MAX_ETA_HEP17};
0180 
0181   const MEbinning eta_binning_hem17_{N_ETA_HEP17, -MAX_ETA_HEP17, MIN_ETA_HEP17};
0182 };
0183 
0184 //
0185 // constants, enums and typedefs
0186 //
0187 
0188 //
0189 // static data member definitions
0190 //
0191 hltPlot HLTObjectsMonitor::getPlotPSet(edm::ParameterSet pset) {
0192   return hltPlot{std::make_pair<MonitorElement*, bool>(nullptr, false),
0193                  std::make_pair<MonitorElement*, bool>(nullptr, pset.getParameter<bool>("displayInPrimary_eta")),
0194                  std::make_pair<MonitorElement*, bool>(nullptr, pset.getParameter<bool>("displayInPrimary_phi")),
0195                  std::make_pair<MonitorElement*, bool>(nullptr, pset.getParameter<bool>("displayInPrimary_pt")),
0196                  std::make_pair<MonitorElement*, bool>(nullptr, pset.getParameter<bool>("displayInPrimary_mass")),
0197                  std::make_pair<MonitorElement*, bool>(nullptr, pset.getParameter<bool>("displayInPrimary_energy")),
0198                  std::make_pair<MonitorElement*, bool>(nullptr, pset.getParameter<bool>("displayInPrimary_csv")),
0199                  std::make_pair<MonitorElement*, bool>(nullptr, pset.getParameter<bool>("displayInPrimary_etaVSphi")),
0200                  std::make_pair<MonitorElement*, bool>(nullptr, pset.getParameter<bool>("displayInPrimary_pt_HEP17")),
0201                  std::make_pair<MonitorElement*, bool>(nullptr, pset.getParameter<bool>("displayInPrimary_pt_HEM17")),
0202                  std::make_pair<MonitorElement*, bool>(nullptr, pset.getParameter<bool>("displayInPrimary_MR")),
0203                  std::make_pair<MonitorElement*, bool>(nullptr, pset.getParameter<bool>("displayInPrimary_RSQ")),
0204                  std::make_pair<MonitorElement*, bool>(nullptr, pset.getParameter<bool>("displayInPrimary_dxy")),
0205                  std::make_pair<MonitorElement*, bool>(nullptr, pset.getParameter<bool>("displayInPrimary_dz")),
0206                  std::make_pair<MonitorElement*, bool>(nullptr, pset.getParameter<bool>("displayInPrimary_dimass")),
0207                  std::make_pair<MonitorElement*, bool>(nullptr, false),
0208                  std::make_pair<MonitorElement*, bool>(nullptr, false),
0209                  std::make_pair<MonitorElement*, bool>(nullptr, false),
0210                  pset.getParameter<std::string>("label"),
0211                  pset.getParameter<std::string>("pathNAME"),
0212                  -1,
0213                  pset.getParameter<std::string>("moduleNAME"),
0214                  pset.getParameter<std::string>("xTITLE"),
0215                  pset.getParameter<std::vector<double>>("etaBINNING"),
0216                  pset.getParameter<std::vector<double>>("ptBINNING"),
0217                  pset.getParameter<std::vector<double>>("phiBINNING"),
0218                  pset.getParameter<std::vector<double>>("massBINNING"),
0219                  pset.getParameter<std::vector<double>>("dxyBINNING"),
0220                  pset.getParameter<std::vector<double>>("dzBINNING"),
0221                  pset.getParameter<std::vector<double>>("dimassBINNING"),
0222                  pset.getUntrackedParameter<bool>("doPlot2D", false),
0223                  pset.getUntrackedParameter<bool>("doPlotETA", true),
0224                  pset.getUntrackedParameter<bool>("doPlotMASS", false),
0225                  pset.getUntrackedParameter<bool>("doPlotENERGY", false),
0226                  pset.getUntrackedParameter<bool>("doPlotHEP17", true),
0227                  pset.getUntrackedParameter<bool>("doPlotCSV", false),
0228                  pset.getUntrackedParameter<bool>("doCALO", false),
0229                  pset.getUntrackedParameter<bool>("doPF", false),
0230                  pset.getUntrackedParameter<bool>("doPlotRazor", false),
0231                  pset.getUntrackedParameter<bool>("doPlotDXY", false),
0232                  pset.getUntrackedParameter<bool>("doPlotDZ", false),
0233                  pset.getUntrackedParameter<bool>("doPlotDiMass", false)};
0234 }
0235 
0236 void HLTObjectsMonitor::getPSet() {
0237   for (const auto& pset : plotPSETS_)
0238     hltPlots_.push_back(getPlotPSet(pset));
0239 }
0240 
0241 bool HLTObjectsMonitor::isHEP17(double eta, double phi) {
0242   if ((eta >= MIN_ETA_HEP17 && eta <= MAX_ETA_HEP17) && (phi >= MIN_PHI_HEP17 && phi <= MAX_PHI_HEP17))
0243     return true;
0244   else
0245     return false;
0246 }
0247 bool HLTObjectsMonitor::isHEM17(double eta, double phi) {
0248   if ((eta >= -MAX_ETA_HEP17 && eta <= -MIN_ETA_HEP17) && (phi >= MIN_PHI_HEP17 && phi <= MAX_PHI_HEP17))
0249     return true;
0250   else
0251     return false;
0252 }
0253 //
0254 // constructors and destructor
0255 //
0256 HLTObjectsMonitor::HLTObjectsMonitor(const edm::ParameterSet& iConfig)
0257     : TopFolder_(iConfig.getParameter<std::string>("TopFolder")),
0258       label_(iConfig.getParameter<std::string>("label")),
0259       processName_(iConfig.getParameter<std::string>("processName")),
0260       plotPSETS_(iConfig.getParameter<std::vector<edm::ParameterSet>>("plots")),
0261       debug_(iConfig.getUntrackedParameter<bool>("debug", false)),
0262       triggerResultsToken_(consumes<edm::TriggerResults>(iConfig.getParameter<edm::InputTag>("TriggerResults"))),
0263       triggerEventToken_(consumes<trigger::TriggerEvent>(iConfig.getParameter<edm::InputTag>("TriggerSummary"))),
0264       beamSpot_(iConfig.getParameter<edm::InputTag>("beamspot")),
0265       beamSpotToken_(consumes<reco::BeamSpot>(beamSpot_)),
0266       caloJetBTagsToken_(consumes<reco::JetTagCollection>(iConfig.getParameter<edm::InputTag>("caloJetBTags"))),
0267       pfJetBTagsToken_(consumes<reco::JetTagCollection>(iConfig.getParameter<edm::InputTag>("pfJetBTags"))),
0268       muCandidates_(iConfig.getParameter<edm::InputTag>("muCandidates")),
0269       muCandidatesToken_(consumes<std::vector<reco::RecoChargedCandidate>>(muCandidates_)),
0270       eleCandidates_(iConfig.getParameter<edm::InputTag>("eleCandidates")),
0271       eleCandidatesToken_(consumes<std::vector<reco::RecoChargedCandidate>>(eleCandidates_)) {
0272   getPSet();
0273 
0274   //now do what ever initialization is needed
0275   mainFolder_ = TopFolder_ + "/MainShifter";
0276   backupFolder_ = TopFolder_ + "/Backup";
0277 
0278   //set Token(s)
0279 
0280   double step = 2 * MAX_PHI / double(N_PHI);
0281   for (int i = 0; i <= N_PHI; i++)
0282     phi_variable_binning_.push_back(-MAX_PHI + step * i);
0283 }
0284 
0285 //
0286 // member functions
0287 //
0288 
0289 // ------------ method called for each event  ------------
0290 void HLTObjectsMonitor::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0291   //  if ( debug_ )
0292   //    std::cout << "[HLTObjectsMonitor::analyze]" << std::endl;
0293 
0294   // access trigger results
0295   edm::Handle<edm::TriggerResults> triggerResults;
0296   iEvent.getByToken(triggerResultsToken_, triggerResults);
0297   if (!triggerResults.isValid())
0298     return;
0299 
0300   edm::Handle<trigger::TriggerEvent> triggerEvent;
0301   iEvent.getByToken(triggerEventToken_, triggerEvent);
0302   if (!triggerEvent.isValid())
0303     return;
0304 
0305   edm::Handle<reco::JetTagCollection> caloJetBTags;
0306   iEvent.getByToken(caloJetBTagsToken_, caloJetBTags);
0307 
0308   edm::Handle<reco::JetTagCollection> pfJetBTags;
0309   iEvent.getByToken(pfJetBTagsToken_, pfJetBTags);
0310 
0311   edm::Handle<std::vector<reco::RecoChargedCandidate>> muCandidates;
0312   iEvent.getByToken(muCandidatesToken_, muCandidates);
0313 
0314   edm::Handle<std::vector<reco::RecoChargedCandidate>> eleCandidates;
0315   iEvent.getByToken(eleCandidatesToken_, eleCandidates);
0316 
0317   edm::Handle<reco::BeamSpot> beamspot;
0318   iEvent.getByToken(beamSpotToken_, beamspot);
0319 
0320   // loop over path
0321   int ibin = -1;
0322   std::vector<bool> plottedPathIndices(plotNamesToBins_.size(), false);
0323   for (auto& plot : hltPlots_) {
0324     ibin++;
0325     if (plot.pathIDX <= 0)
0326       continue;
0327 
0328     if (triggerResults->accept(plot.pathIDX)) {
0329       //We only want to fill this once per pathNAME per event
0330       auto index = plotNamesToBins_[plot.pathNAME];
0331       if (not plottedPathIndices[index]) {
0332         plottedPathIndices[index] = true;
0333         if (debug_)
0334           std::cout << plot.pathNAME << " --> bin: " << ibin << std::endl;
0335         eventsPlot_->Fill(index);
0336       }
0337       const trigger::TriggerObjectCollection objects = triggerEvent->getObjects();
0338       if (hltConfig_.saveTags(plot.moduleNAME)) {
0339         if (debug_)
0340           std::cout << "objects: " << objects.size() << std::endl;
0341 
0342         bool moduleFOUND = false;
0343         std::vector<std::string> moduleNames = hltConfig_.moduleLabels(plot.pathIDX);
0344         for (const auto& module : moduleNames) {
0345           if (module == plot.moduleNAME)
0346             moduleFOUND = true;
0347         }
0348         if (debug_)
0349           std::cout << plot.moduleNAME << (moduleFOUND ? "" : "NOT") << " found in the list of modules" << std::endl;
0350 
0351         if (debug_)
0352           for (const auto& module : moduleNames) {
0353             unsigned int idx = triggerEvent->filterIndex(edm::InputTag(module, "", processName_));
0354             std::cout << "module: " << module;
0355             if (idx < triggerEvent->sizeFilters())
0356               std::cout << " --> " << idx;
0357             std::cout << std::endl;
0358           }
0359         //
0360         // trigger accepted and collection w/ objects is available
0361         edm::InputTag moduleName = edm::InputTag(plot.moduleNAME, "", processName_);
0362         unsigned int moduleIDX = triggerEvent->filterIndex(moduleName);
0363         if (debug_)
0364           std::cout << "moduleNAME: " << plot.moduleNAME << " --> " << moduleIDX << std::endl;
0365 
0366         if (moduleIDX >= triggerEvent->sizeFilters()) {
0367           LogDebug("HLTObjectsMonitor")
0368               << plot.pathNAME << " " << plot.moduleNAME
0369               << " is not available ! please, fix update DQM/HLTEvF/python/HLTObjectsMonitor_cfi.py";
0370           return;
0371         }
0372 
0373         const trigger::Keys& keys = triggerEvent->filterKeys(moduleIDX);
0374         if (debug_)
0375           std::cout << "keys: " << keys.size() << std::endl;
0376 
0377         plot.nME.first->Fill(keys.size());
0378 
0379         double MR = 0.;
0380         double RSQ = 0.;
0381         for (const auto& key : keys) {
0382           double pt = objects[key].pt();
0383           double eta = objects[key].eta();
0384           double phi = objects[key].phi();
0385           double mass = objects[key].mass();
0386           double energy = objects[key].energy();
0387           int id = objects[key].id();
0388           if (debug_)
0389             std::cout << "object ID " << id << " mass: " << mass << std::endl;
0390 
0391           // single-object plots
0392           plot.ptME.first->Fill(pt);
0393           if (plot.doPlotETA)
0394             plot.etaME.first->Fill(eta);
0395           plot.phiME.first->Fill(phi);
0396 
0397           if (plot.doPlotCSV) {
0398             if (plot.doCALO) {
0399               if (!caloJetBTags.isValid())
0400                 plot.csvME.first->Fill(-10.);
0401               else {
0402                 for (auto it = caloJetBTags->begin(); it != caloJetBTags->end(); ++it) {
0403                   double dR = deltaR(eta, phi, it->first->eta(), it->first->phi());
0404                   if (debug_)
0405                     std::cout << "[HLTObjectsMonitor::analyze] deltaR: " << dR << " matched ? "
0406                               << (dR <= 0.4 ? "YEAP" : "NOPE") << std::endl;
0407                   plot.csvME.first->Fill(it->second);
0408                 }
0409               }
0410 
0411             } else if (plot.doPF) {
0412               if (!pfJetBTags.isValid())
0413                 plot.csvME.first->Fill(-10.);
0414               else {
0415                 for (auto it = pfJetBTags->begin(); it != pfJetBTags->end(); ++it) {
0416                   double dR = deltaR(eta, phi, it->first->eta(), it->first->phi());
0417                   if (debug_)
0418                     std::cout << "[HLTObjectsMonitor::analyze] deltaR: " << dR << " matched ? "
0419                               << (dR <= 0.4 ? "YEAP" : "NOPE") << std::endl;
0420                   plot.csvME.first->Fill(it->second);
0421                 }
0422               }
0423             }
0424           }
0425           if (plot.doPlotMASS)
0426             plot.massME.first->Fill(mass);
0427           if (plot.doPlotENERGY)
0428             plot.energyME.first->Fill(energy);
0429           if (plot.doPlot2D)
0430             plot.etaVSphiME.first->Fill(eta, phi);
0431           if (plot.doPlotHEP17) {
0432             if (isHEP17(eta, phi))
0433               plot.ptMEhep17.first->Fill(pt);
0434             if (isHEM17(eta, phi))
0435               plot.ptMEhem17.first->Fill(pt);
0436           }
0437 
0438           if (id == 0) {             //the MET object containing MR and Rsq will show up with ID = 0
0439             MR = objects[key].px();  //razor variables stored in dummy reco::MET objects
0440             RSQ = objects[key].py();
0441           }
0442 
0443           if (plot.doPlotDXY) {
0444             double dxy = -99.;
0445             if (abs(id) == 13)
0446               dxy = dxyFinder(eta, phi, muCandidates, beamspot, 0.1);  // dRcut = 0.1
0447             else
0448               dxy = dxyFinder(eta, phi, eleCandidates, beamspot, 0.1);  // dRcut = 0.1
0449             plot.dxyME.first->Fill(dxy);
0450           }
0451         }  // end loop on keys
0452         if (plot.doPlotRazor) {
0453           plot.mrME.first->Fill(MR);
0454           plot.rsqME.first->Fill(RSQ);
0455         }
0456 
0457         if (keys.size() < 2) {
0458           if (plot.doPlotDiMass || plot.doPlotDZ)
0459             LogDebug("HLTObjectsMonitor") << plot.pathNAME << " " << plot.moduleNAME << " # objects is (" << keys.size()
0460                                           << ") less than 2 ! you probably want to either change the moduleNAME or "
0461                                              "switch off di-object system plots (doPlotDZ: "
0462                                           << plot.doPlotDZ << " doPlotDiMass: " << plot.doPlotDiMass
0463                                           << ") in DQM/HLTEvF/python/HLTObjectsMonitor_cfi.py)";
0464         } else {
0465           for (const auto& key : keys) {
0466             double pt = objects[key].pt();
0467             double eta = objects[key].eta();
0468             double phi = objects[key].phi();
0469             int id = objects[key].id();
0470 
0471             unsigned int kCnt0 = 0;
0472 
0473             TLorentzVector v1;
0474             if (abs(id) == 13)  // check if it is a muon
0475               v1.SetPtEtaPhiM(pt, eta, phi, MASS_MU);
0476             else
0477               v1.SetPtEtaPhiM(pt, eta, phi, 0);
0478 
0479             unsigned int kCnt1 = 0;
0480             for (const auto& key1 : keys) {
0481               if (key != key1 &&
0482                   kCnt1 > kCnt0) {  // avoid filling hists with same objs && avoid double counting separate objs
0483 
0484                 double pt2 = objects[key1].pt();
0485                 double eta2 = objects[key1].eta();
0486                 double phi2 = objects[key1].phi();
0487                 int id2 = objects[key1].id();
0488 
0489                 double dR = deltaR(eta, phi, eta2, phi2);
0490                 plot.dRME.first->Fill(dR);
0491                 plot.dRetaVSphiME.first->Fill(eta, phi, dR);
0492 
0493                 int q1 = (id == 0 ? 0 : id / abs(id));
0494                 int q2 = (id2 == 0 ? 0 : id2 / abs(id2));
0495                 int q1q2 = q1 * q2;
0496                 plot.q1q2ME.first->Fill(q1q2);
0497 
0498                 if ((id + id2) == 0) {  // check di-object system charge and flavor
0499 
0500                   TLorentzVector v2;
0501                   if (abs(id2) == 13)  // check if it is a muon
0502                     v2.SetPtEtaPhiM(pt2, eta2, phi2, MASS_MU);
0503                   else
0504                     v2.SetPtEtaPhiM(pt2, eta2, phi2, 0);
0505 
0506                   if (plot.doPlotDiMass) {
0507                     TLorentzVector v = v1 + v2;
0508                     plot.dimassME.first->Fill(v.M());
0509                   }
0510 
0511                   if (plot.doPlotDZ) {
0512                     double dz = -99.;
0513                     if (abs(id) == 13)
0514                       dz = dzFinder(eta, phi, eta2, phi2, muCandidates, 0.1);  // dRcut = 0.1
0515                     else
0516                       dz = dzFinder(eta, phi, eta2, phi2, eleCandidates, 0.1);  // dRcut = 0.1
0517                     plot.dzME.first->Fill(dz);
0518                   }
0519                 }
0520               }
0521               kCnt1++;
0522             }
0523             kCnt0++;
0524           }
0525         }
0526       }
0527     }
0528   }
0529 }
0530 
0531 // ------------ method called when starting to processes a run  ------------
0532 void HLTObjectsMonitor::dqmBeginRun(edm::Run const& iRun, edm::EventSetup const& iSetup) {
0533   bool changed = true;
0534   if (hltConfig_.init(iRun, iSetup, processName_, changed))
0535     if (debug_)
0536       std::cout << "[HLTObjectsMonitor::dqmBeginRun] extracting HLTconfig" << std::endl;
0537 
0538   //get path indicies from menu
0539   std::string pathName_noVersion;
0540   std::vector<std::string> triggerPaths = hltConfig_.triggerNames();
0541 
0542   if (debug_)
0543     std::cout << "[HLTObjectsMonitor::dqmBeginRun] triggerPaths: " << triggerPaths.size() << " <--> "
0544               << hltPlots_.size() << std::endl;
0545 
0546   for (const auto& pathName : triggerPaths) {
0547     if (debug_)
0548       std::cout << "[HLTObjectsMonitor::dqmBeginRun] " << pathName << std::endl;
0549 
0550     pathName_noVersion = hltConfig_.removeVersion(pathName);
0551     //    std::cout << "pathName: " << pathName << " --> " << pathName_noVersion << std::endl;
0552     for (auto& plot : hltPlots_) {
0553       if (plot.pathNAME == pathName_noVersion) {
0554         plot.pathIDX = hltConfig_.triggerIndex(pathName);
0555         // check that the index makes sense, otherwise force pathIDX = -1
0556         if (plot.pathIDX <= 0 || plot.pathIDX == int(triggerPaths.size()))
0557           plot.pathIDX = -1;
0558       }
0559     }
0560   }
0561 
0562   if (debug_) {
0563     for (const auto& plot : hltPlots_)
0564       std::cout << "plot: " << plot.pathNAME << " --> pathIDX: " << plot.pathIDX << std::endl;
0565     std::cout << "[HLTObjectsMonitor::dqmBeginRun] DONE" << std::endl;
0566   }
0567 }
0568 
0569 void HLTObjectsMonitor::bookHistograms(DQMStore::IBooker& ibooker,
0570                                        edm::Run const& iRun,
0571                                        edm::EventSetup const& iSetup) {
0572   if (debug_)
0573     std::cout << "[HLTObjectsMonitor::bookHistograms]" << std::endl;
0574 
0575   ibooker.setCurrentFolder(TopFolder_);
0576 
0577   std::string name = "eventsPerPath_" + label_;
0578   std::string title = " events per path";
0579 
0580   //We must avoid repeating the same pathNAME
0581   {
0582     std::unordered_map<std::string, int> uniqueNames;
0583     for (auto const& p : hltPlots_) {
0584       plotNamesToBins_[p.pathNAME] = -1;
0585     }
0586     int nbins = plotNamesToBins_.size();
0587     eventsPlot_ = ibooker.book1D(name, title, nbins, -0.5, double(nbins) - 0.5);
0588     eventsPlot_->setAxisTitle("HLT path");
0589     int i = 0;
0590     //keep the bin order the same as hltPlots_
0591     for (auto const& p : hltPlots_) {
0592       //only add a bin if this is the first time we've seen the name
0593       if (-1 == plotNamesToBins_[p.pathNAME]) {
0594         plotNamesToBins_[p.pathNAME] = ++i;
0595         eventsPlot_->setBinLabel(i, p.pathNAME);
0596         if (debug_)
0597           std::cout << p.pathNAME << " --> bin: " << i << std::endl;
0598       }
0599     }
0600   }
0601 
0602   for (auto& plot : hltPlots_) {
0603     if (debug_)
0604       std::cout << "booking plots for " << plot.label << std::endl;
0605 
0606     if (plot.pathIDX <= 0) {
0607       LogDebug("HLTObjectsMonitor") << plot.pathNAME
0608                                     << " is not available in the HLT menu ! no plots are going to be booked for it "
0609                                        "(update DQM/HLTEvF/python/HLTObjectsMonitor_cfi.py)";
0610       continue;
0611     }
0612     if (debug_)
0613       std::cout << "booking histograms for " << plot.pathNAME << std::endl;
0614 
0615     {
0616       if (plot.nME.second)
0617         ibooker.setCurrentFolder(mainFolder_);
0618       else
0619         ibooker.setCurrentFolder(backupFolder_);
0620 
0621       name = plot.label + "_nobjects";
0622       title = plot.pathNAME + " # objects";
0623       plot.nME.first = ibooker.book1D(name, title, 20, -0.5, 19.5);
0624       plot.nME.first->setAxisTitle(plot.xTITLE + " # objects");
0625     }
0626 
0627     if (plot.ptME.second)
0628       ibooker.setCurrentFolder(mainFolder_);
0629     else
0630       ibooker.setCurrentFolder(backupFolder_);
0631 
0632     name = plot.label + "_pt";
0633     title = plot.pathNAME + " p_T";
0634     int nbins = (plot.ptBINNING).size() - 1;
0635     std::vector<float> fbinning((plot.ptBINNING).begin(), (plot.ptBINNING).end());
0636     float* arr = &fbinning[0];
0637     plot.ptME.first = ibooker.book1D(name, title, nbins, arr);
0638     plot.ptME.first->setAxisTitle(plot.xTITLE + " p_{T} [GeV]");
0639 
0640     {
0641       if (plot.phiME.second)
0642         ibooker.setCurrentFolder(mainFolder_);
0643       else
0644         ibooker.setCurrentFolder(backupFolder_);
0645 
0646       name = plot.label + "_phi";
0647       title = plot.pathNAME + " #phi";
0648       int nbins = (plot.phiBINNING).size() - 1;
0649       std::vector<float> fbinning((plot.phiBINNING).begin(), (plot.phiBINNING).end());
0650       float* arr = &fbinning[0];
0651       plot.phiME.first = ibooker.book1D(name, title, nbins, arr);
0652       plot.phiME.first->setAxisTitle(plot.xTITLE + " #phi [rad]");
0653     }
0654 
0655     if (plot.doPlotETA) {
0656       if (plot.etaME.second)
0657         ibooker.setCurrentFolder(mainFolder_);
0658       else
0659         ibooker.setCurrentFolder(backupFolder_);
0660 
0661       name = plot.label + "_eta";
0662       title = plot.pathNAME + " #eta";
0663       int nbins = (plot.etaBINNING).size() - 1;
0664       std::vector<float> fbinning((plot.etaBINNING).begin(), (plot.etaBINNING).end());
0665       float* arr = &fbinning[0];
0666 
0667       plot.etaME.first = ibooker.book1D(name, title, nbins, arr);
0668       plot.etaME.first->setAxisTitle(plot.xTITLE + " #eta");
0669     }
0670 
0671     if (plot.doPlotMASS) {
0672       if (plot.massME.second)
0673         ibooker.setCurrentFolder(mainFolder_);
0674       else
0675         ibooker.setCurrentFolder(backupFolder_);
0676 
0677       name = plot.label + "_mass";
0678       title = plot.pathNAME + " mass";
0679       int nbins = (plot.massBINNING).size() - 1;
0680       std::vector<float> fbinning((plot.massBINNING).begin(), (plot.massBINNING).end());
0681       float* arr = &fbinning[0];
0682 
0683       plot.massME.first = ibooker.book1D(name, title, nbins, arr);
0684       plot.massME.first->setAxisTitle(plot.xTITLE + " mass [GeV]");
0685     }
0686 
0687     if (plot.doPlotENERGY) {
0688       if (plot.energyME.second)
0689         ibooker.setCurrentFolder(mainFolder_);
0690       else
0691         ibooker.setCurrentFolder(backupFolder_);
0692 
0693       name = plot.label + "_energy";
0694       title = plot.pathNAME + " energy";
0695       int nbins = (plot.ptBINNING).size() - 1;
0696       std::vector<float> fbinning((plot.ptBINNING).begin(), (plot.ptBINNING).end());
0697       float* arr = &fbinning[0];
0698 
0699       plot.energyME.first = ibooker.book1D(name, title, nbins, arr);
0700       plot.energyME.first->setAxisTitle(plot.xTITLE + " energy [GeV]");
0701     }
0702 
0703     if (plot.doPlotCSV) {
0704       if (plot.csvME.second)
0705         ibooker.setCurrentFolder(mainFolder_);
0706       else
0707         ibooker.setCurrentFolder(backupFolder_);
0708 
0709       name = plot.label + "_csv";
0710       title = plot.pathNAME + " CSV";
0711 
0712       plot.csvME.first = ibooker.book1D(name, title, csv_binning_.nbins, csv_binning_.xmin, csv_binning_.xmax);
0713       plot.csvME.first->setAxisTitle(plot.xTITLE + " CSV discriminator");
0714     }
0715 
0716     if (plot.doPlot2D) {
0717       if (plot.etaVSphiME.second)
0718         ibooker.setCurrentFolder(mainFolder_);
0719       else
0720         ibooker.setCurrentFolder(backupFolder_);
0721 
0722       name = plot.label + "_etaVSphi";
0723       title = plot.pathNAME + " #eta vs #phi";
0724       int nbinsX = (plot.etaBINNING).size() - 1;
0725       std::vector<float> fbinningX((plot.etaBINNING).begin(), (plot.etaBINNING).end());
0726       float* arrX = &fbinningX[0];
0727       int nbinsY = (plot.phiBINNING).size() - 1;
0728       ;
0729       std::vector<float> fbinningY((plot.phiBINNING).begin(), (plot.phiBINNING).end());
0730       float* arrY = &fbinningY[0];
0731       plot.etaVSphiME.first = ibooker.book2D(name, title, nbinsX, arrX, nbinsY, arrY);
0732       plot.etaVSphiME.first->setAxisTitle(plot.xTITLE + " #eta", 1);
0733       plot.etaVSphiME.first->setAxisTitle(plot.xTITLE + " #phi", 2);
0734     }
0735 
0736     if (plot.doPlotHEP17) {
0737       if (plot.ptMEhep17.second)
0738         ibooker.setCurrentFolder(mainFolder_);
0739       else
0740         ibooker.setCurrentFolder(backupFolder_);
0741 
0742       int nbins = (plot.ptBINNING).size() - 1;
0743       std::vector<float> fbinning((plot.ptBINNING).begin(), (plot.ptBINNING).end());
0744       float* arr = &fbinning[0];
0745 
0746       name = plot.label + "_pt_HEP17";
0747       title = plot.pathNAME + " p_{T} HEP17";
0748       plot.ptMEhep17.first = ibooker.book1D(name, title, nbins, arr);
0749       plot.ptMEhep17.first->setAxisTitle(plot.xTITLE + " p_{T} [GeV]", 1);
0750 
0751       if (plot.ptMEhem17.second)
0752         ibooker.setCurrentFolder(mainFolder_);
0753       else
0754         ibooker.setCurrentFolder(backupFolder_);
0755 
0756       name = plot.label + "_pt_HEM17";
0757       title = plot.pathNAME + " p_{T} HEM17";
0758       plot.ptMEhem17.first = ibooker.book1D(name, title, nbins, arr);
0759       plot.ptMEhem17.first->setAxisTitle(plot.xTITLE + " p_{T} [GeV]", 1);
0760     }
0761 
0762     if (plot.doPlotRazor) {
0763       if (plot.mrME.second)
0764         ibooker.setCurrentFolder(mainFolder_);
0765       else
0766         ibooker.setCurrentFolder(backupFolder_);
0767 
0768       name = plot.label + "_mr";
0769       title = plot.pathNAME + " M_{R}";
0770       plot.mrME.first = ibooker.book1D(name, title, 100, 0., 100.);
0771       plot.mrME.first->setAxisTitle(plot.xTITLE + " M_{R} [GeV]", 1);
0772 
0773       if (plot.rsqME.second)
0774         ibooker.setCurrentFolder(mainFolder_);
0775       else
0776         ibooker.setCurrentFolder(backupFolder_);
0777 
0778       name = plot.label + "_rsq";
0779       title = plot.pathNAME + " RSQ";
0780       plot.rsqME.first = ibooker.book1D(name, title, 100, 0., 100.);
0781       plot.rsqME.first->setAxisTitle(plot.xTITLE + " RSQ [GeV]", 1);
0782     }
0783 
0784     if (plot.doPlotDXY) {
0785       if (plot.dxyME.second)
0786         ibooker.setCurrentFolder(mainFolder_);
0787       else
0788         ibooker.setCurrentFolder(backupFolder_);
0789 
0790       name = plot.label + "_dxy";
0791       title = plot.pathNAME + " d_{xy}";
0792       int nbins = (plot.dxyBINNING).size() - 1;
0793       std::vector<float> fbinning((plot.dxyBINNING).begin(), (plot.dxyBINNING).end());
0794       float* arr = &fbinning[0];
0795       plot.dxyME.first = ibooker.book1D(name, title, nbins, arr);
0796       plot.dxyME.first->setAxisTitle(plot.xTITLE + " d_{xy} [cm]");
0797     }
0798 
0799     if (plot.doPlotDZ) {
0800       if (plot.dzME.second)
0801         ibooker.setCurrentFolder(mainFolder_);
0802       else
0803         ibooker.setCurrentFolder(backupFolder_);
0804 
0805       name = plot.label + "_dz";
0806       title = plot.pathNAME + " d_{z}";
0807       int nbins = (plot.dzBINNING).size() - 1;
0808       std::vector<float> fbinning((plot.dzBINNING).begin(), (plot.dzBINNING).end());
0809       float* arr = &fbinning[0];
0810       plot.dzME.first = ibooker.book1D(name, title, nbins, arr);
0811       plot.dzME.first->setAxisTitle(plot.xTITLE + " d_{z} [cm]");
0812     }
0813 
0814     if (plot.dRME.second)
0815       ibooker.setCurrentFolder(mainFolder_);
0816     else
0817       ibooker.setCurrentFolder(backupFolder_);
0818 
0819     name = plot.label + "_dR";
0820     title = plot.pathNAME + " di-object dR";
0821     plot.dRME.first = ibooker.book1D(name, title, 50, 0., 5.);
0822     plot.dRME.first->setAxisTitle(plot.xTITLE + " dR_{obj,obj}");
0823 
0824     if (plot.dRetaVSphiME.second)
0825       ibooker.setCurrentFolder(mainFolder_);
0826     else
0827       ibooker.setCurrentFolder(backupFolder_);
0828 
0829     name = plot.label + "_dR_etaVSphi";
0830     title = plot.pathNAME + " di-object dR in the #eta-#phi plane (of 1st obj)";
0831     plot.dRetaVSphiME.first = ibooker.bookProfile2D(name, title, 60, -3., 3., 64, -3.2, 3.2, 0., 5.);
0832     plot.dRetaVSphiME.first->setAxisTitle(plot.xTITLE + " #eta", 1);
0833     plot.dRetaVSphiME.first->setAxisTitle(plot.xTITLE + " #phi", 2);
0834     plot.dRetaVSphiME.first->setAxisTitle(plot.xTITLE + " dR_{obj,obj}", 3);
0835 
0836     if (plot.q1q2ME.second)
0837       ibooker.setCurrentFolder(mainFolder_);
0838     else
0839       ibooker.setCurrentFolder(backupFolder_);
0840 
0841     name = plot.label + "_q1q2";
0842     title = plot.pathNAME + " di-object q1xq2";
0843     plot.q1q2ME.first = ibooker.book1D(name, title, 3, -1., 1.);
0844     plot.q1q2ME.first->setAxisTitle(plot.xTITLE + " q_{obj1} x q_{obj2}");
0845 
0846     if (plot.doPlotDiMass) {
0847       if (plot.dimassME.second)
0848         ibooker.setCurrentFolder(mainFolder_);
0849       else
0850         ibooker.setCurrentFolder(backupFolder_);
0851 
0852       name = plot.label + "_dimass";
0853       title = plot.pathNAME + " di-object mass";
0854       int nbins = (plot.dimassBINNING).size() - 1;
0855       std::vector<float> fbinning((plot.dimassBINNING).begin(), (plot.dimassBINNING).end());
0856       float* arr = &fbinning[0];
0857       plot.dimassME.first = ibooker.book1D(name, title, nbins, arr);
0858       plot.dimassME.first->setAxisTitle(plot.xTITLE + " m_{obj,obj} [GeV]");
0859     }
0860   }
0861 
0862   if (debug_)
0863     std::cout << "[HLTObjectsMonitor::bookHistograms] DONE" << std::endl;
0864 }
0865 
0866 double HLTObjectsMonitor::dxyFinder(double eta,
0867                                     double phi,
0868                                     edm::Handle<reco::RecoChargedCandidateCollection> candidates,
0869                                     edm::Handle<reco::BeamSpot> beamspot,
0870                                     double dRcut = 0.1) {
0871   double dxy = -99.;
0872   if (!candidates.isValid()) {
0873     LogDebug("HLTObjectsMonitor") << "either " << muCandidates_ << " or " << eleCandidates_
0874                                   << " is not valid ! please, update DQM/HLTEvF/python/HLTObjectsMonitor_cfi.py"
0875                                   << " by switching OFF doPlotDXY or updating the InputTag collection";
0876     return dxy;
0877   }
0878   if (!beamspot.isValid()) {
0879     LogDebug("HLTObjectsMonitor") << beamSpot_
0880                                   << " is not valid ! please, update DQM/HLTEvF/python/HLTObjectsMonitor_cfi.py"
0881                                   << " by switching OFF doPlotDXY or updating the InputTag collection";
0882     return dxy;
0883   }
0884 
0885   bool matched = false;
0886   for (reco::RecoChargedCandidateCollection::const_iterator candidate = candidates->begin();
0887        candidate != candidates->end();
0888        ++candidate) {
0889     if (deltaR(eta, phi, candidate->eta(), candidate->phi()) < dRcut) {
0890       matched = true;
0891       dxy = (-(candidate->vx() - beamspot->x0()) * candidate->py() +
0892              (candidate->vy() - beamspot->y0()) * candidate->px()) /
0893             candidate->pt();
0894       break;
0895     }
0896   }
0897   if (!matched)
0898     edm::LogWarning("HLTObjectsMonitor") << "trigger object does not match ( dR > " << dRcut
0899                                          << ") to any of the candidates in either " << muCandidates_ << " or "
0900                                          << eleCandidates_;
0901 
0902   return dxy;
0903 }
0904 
0905 double HLTObjectsMonitor::dzFinder(double eta1,
0906                                    double phi1,
0907                                    double eta2,
0908                                    double phi2,
0909                                    edm::Handle<reco::RecoChargedCandidateCollection> candidates,
0910                                    double dRcut = 0.1) {
0911   double dz = -99.;
0912   if (!candidates.isValid()) {
0913     LogDebug("HLTObjectsMonitor") << "either " << muCandidates_ << " or " << eleCandidates_
0914                                   << " is not valid ! please, update DQM/HLTEvF/python/HLTObjectsMonitor_cfi.py"
0915                                   << " by switching OFF doPlotDZ or updating the InputTag collection";
0916     return dz;
0917   }
0918 
0919   const reco::RecoChargedCandidate* cand1;
0920   const reco::RecoChargedCandidate* cand2;
0921   bool matched1 = false;
0922   bool matched2 = false;
0923   for (reco::RecoChargedCandidateCollection::const_iterator candidate = candidates->begin();
0924        candidate != candidates->end();
0925        ++candidate) {
0926     if (deltaR(eta1, phi1, candidate->eta(), candidate->phi()) < dRcut) {
0927       matched1 = true;
0928       cand1 = &*candidate;
0929       if (debug_)
0930         std::cout << "cand1: " << cand1->pt() << " " << cand1->eta() << " " << cand1->phi() << std::endl;
0931       break;
0932     }
0933   }
0934   if (!matched1) {
0935     LogDebug("HLTObjectsMonitor") << "trigger object1 does not match ( dR > " << dRcut
0936                                   << ") to any of the candidates in either " << muCandidates_ << " or "
0937                                   << eleCandidates_;
0938     return dz;
0939   }
0940 
0941   for (reco::RecoChargedCandidateCollection::const_iterator candidate = candidates->begin();
0942        candidate != candidates->end();
0943        ++candidate) {
0944     if (debug_) {
0945       std::cout << "candidate: " << candidate->pt() << " cand1: " << cand1->pt() << std::endl;
0946       std::cout << "candidate: " << candidate->eta() << " cand1: " << cand1->eta() << std::endl;
0947       std::cout << "candidate: " << candidate->phi() << " cand1: " << cand1->phi() << std::endl;
0948     }
0949     if (&*candidate == cand1)
0950       continue;
0951 
0952     if (deltaR(eta2, phi2, candidate->eta(), candidate->phi()) < dRcut) {
0953       matched2 = true;
0954       cand2 = &*candidate;
0955       if (debug_)
0956         std::cout << "cand2: " << cand2->pt() << " " << cand2->eta() << " " << cand2->phi() << std::endl;
0957       break;
0958     }
0959   }
0960   if (!matched2) {
0961     LogDebug("HLTObjectsMonitor") << "trigger object2 does not match ( dR > " << dRcut
0962                                   << ") to any of the candidates in either " << muCandidates_ << " or "
0963                                   << eleCandidates_;
0964     return dz;
0965   }
0966 
0967   dz = cand1->vz() - cand2->vz();
0968   return dz;
0969 }
0970 
0971 //define this as a plug-in
0972 DEFINE_FWK_MODULE(HLTObjectsMonitor);