Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-11 04:32:52

0001 #include <string>
0002 #include <vector>
0003 
0004 #include "FWCore/Utilities/interface/EDGetToken.h"
0005 #include "FWCore/Utilities/interface/ESGetToken.h"
0006 #include "FWCore/Framework/interface/Frameworkfwd.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0009 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0010 #include "FWCore/ParameterSet/interface/Registry.h"
0011 #include "FWCore/ServiceRegistry/interface/Service.h"
0012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0013 #include "FWCore/Framework/interface/MakerMacros.h"
0014 #include "DQMServices/Core/interface/DQMStore.h"
0015 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0016 #include "DQMOffline/Trigger/plugins/TriggerDQMBase.h"
0017 #include "CommonTools/TriggerUtils/interface/GenericTriggerEventFlag.h"
0018 #include "CommonTools/TriggerUtils/interface/PrescaleWeightProvider.h"
0019 #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
0020 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0021 #include "DataFormats/Math/interface/deltaR.h"
0022 #include "DataFormats/Math/interface/deltaPhi.h"
0023 #include "MagneticField/Engine/interface/MagneticField.h"
0024 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0025 #include "TrackingTools/PatternTools/interface/ClosestApproachInRPhi.h"
0026 #include "DataFormats/EgammaCandidates/interface/Photon.h"
0027 #include "DataFormats/EgammaCandidates/interface/PhotonFwd.h"
0028 #include "DataFormats/VertexReco/interface/Vertex.h"
0029 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0030 #include "DataFormats/MuonReco/interface/Muon.h"
0031 #include "DataFormats/MuonReco/interface/MuonFwd.h"
0032 #include "RecoVertex/VertexPrimitives/interface/TransientVertex.h"
0033 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexFitter.h"
0034 #include "DataFormats/TrackReco/interface/Track.h"
0035 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0036 #include "DataFormats/TrackReco/interface/TrackExtra.h"
0037 #include "DataFormats/TrackReco/interface/TrackExtraFwd.h"
0038 #include "DataFormats/HLTReco/interface/TriggerEvent.h"
0039 #include "DataFormats/HLTReco/interface/TriggerEventWithRefs.h"
0040 #include "DataFormats/Common/interface/TriggerResults.h"
0041 #include "DataFormats/HLTReco/interface/TriggerObject.h"
0042 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0043 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
0044 #include "HLTrigger/HLTcore/interface/HLTPrescaleProvider.h"
0045 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
0046 
0047 class BPHMonitor : public DQMEDAnalyzer, public TriggerDQMBase {
0048 public:
0049   typedef dqm::reco::MonitorElement MonitorElement;
0050   typedef dqm::reco::DQMStore DQMStore;
0051 
0052   BPHMonitor(const edm::ParameterSet&);
0053   ~BPHMonitor() throw() override;
0054   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0055 
0056 protected:
0057   void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
0058   void analyze(edm::Event const& iEvent, edm::EventSetup const& iSetup) override;
0059 
0060   template <typename T>
0061   bool matchToTrigger(const std::string& theTriggerName, T t);
0062 
0063   double Prescale(const std::string num,
0064                   const std::string den,
0065                   edm::Event const& iEvent,
0066                   edm::EventSetup const& iSetup,
0067                   HLTPrescaleProvider* hltPrescale_);
0068 
0069 private:
0070   const std::string folderName_;
0071 
0072   const bool requireValidHLTPaths_;
0073   bool hltPathsAreValid_;
0074 
0075   edm::InputTag muoInputTag_;
0076   edm::InputTag bsInputTag_;
0077   edm::InputTag trInputTag_;
0078   edm::InputTag phInputTag_;
0079   edm::InputTag vtxInputTag_;
0080 
0081   edm::EDGetTokenT<reco::MuonCollection> muoToken_;
0082   edm::EDGetTokenT<reco::BeamSpot> bsToken_;
0083   edm::EDGetTokenT<reco::TrackCollection> trToken_;
0084   edm::EDGetTokenT<reco::PhotonCollection> phToken_;
0085   edm::EDGetTokenT<reco::VertexCollection> vtxToken_;
0086 
0087   edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> magneticFieldToken_;
0088 
0089   std::vector<double> pt_variable_binning_;
0090   std::vector<double> dMu_pt_variable_binning_;
0091   std::vector<double> prob_variable_binning_;
0092   MEbinning phi_binning_;
0093   MEbinning eta_binning_;
0094   MEbinning d0_binning_;
0095   MEbinning z0_binning_;
0096   MEbinning dR_binning_;
0097   MEbinning mass_binning_;
0098   MEbinning Bmass_binning_;
0099   MEbinning dca_binning_;
0100   MEbinning ds_binning_;
0101   MEbinning cos_binning_;
0102 
0103   ObjME muPhi_;
0104   ObjME muEta_;
0105   ObjME muPt_;
0106   ObjME mud0_;
0107   ObjME muz0_;
0108 
0109   ObjME mu1Phi_;
0110   ObjME mu1Eta_;
0111   ObjME mu1Pt_;
0112   ObjME mu1d0_;
0113   ObjME mu1z0_;
0114   ObjME mu2Phi_;
0115   ObjME mu2Eta_;
0116   ObjME mu2Pt_;
0117   ObjME mu2d0_;
0118   ObjME mu2z0_;
0119   ObjME mu3Phi_;
0120   ObjME mu3Eta_;
0121   ObjME mu3Pt_;
0122   ObjME mu3d0_;
0123   ObjME mu3z0_;
0124 
0125   ObjME phPhi_;
0126   ObjME phEta_;
0127   ObjME phPt_;
0128   ObjME DiMuPhi_;
0129   ObjME DiMuEta_;
0130   ObjME DiMuPt_;
0131   ObjME DiMuPVcos_;
0132   ObjME DiMuProb_;
0133   ObjME DiMuDS_;
0134   ObjME DiMuDCA_;
0135   ObjME DiMuMass_;
0136   ObjME BMass_;
0137   ObjME DiMudR_;
0138 
0139   std::unique_ptr<GenericTriggerEventFlag> num_genTriggerEventFlag_;
0140   std::unique_ptr<GenericTriggerEventFlag> den_genTriggerEventFlag_;
0141 
0142   HLTPrescaleProvider* hltPrescale_;
0143 
0144   StringCutObjectSelector<reco::Muon, true> muoSelection_;
0145   StringCutObjectSelector<reco::Muon, true> muoSelection_ref;
0146   StringCutObjectSelector<reco::Muon, true> muoSelection_tag;
0147   StringCutObjectSelector<reco::Muon, true> muoSelection_probe;
0148 
0149   int nmuons_;
0150   bool tnp_;
0151   int L3_;
0152   int ptCut_;
0153   int displaced_;
0154   int trOrMu_;
0155   int Jpsi_;
0156   int Upsilon_;
0157   int enum_;
0158   int seagull_;
0159   double maxmass_;
0160   double minmass_;
0161   double maxmassJpsi;
0162   double minmassJpsi;
0163   double maxmassUpsilon;
0164   double minmassUpsilon;
0165   double maxmassTkTk;
0166   double minmassTkTk;
0167   double maxmassJpsiTk;
0168   double minmassJpsiTk;
0169   double kaon_mass;
0170   double mu_mass;
0171   double min_dR;
0172   double max_dR;
0173 
0174   double minprob;
0175   double mincos;
0176   double minDS;
0177   edm::EDGetTokenT<edm::TriggerResults> hltTrigResTag_;
0178   edm::InputTag hltInputTag_1;
0179   edm::EDGetTokenT<trigger::TriggerEvent> hltInputTag_;
0180   std::vector<std::string> hltpaths_num;
0181   std::vector<std::string> hltpaths_den;
0182   StringCutObjectSelector<reco::Track, true> trSelection_;
0183   StringCutObjectSelector<reco::Track, true> trSelection_ref;
0184   StringCutObjectSelector<reco::Candidate::LorentzVector, true> DMSelection_ref;
0185 
0186   edm::Handle<trigger::TriggerEvent> handleTriggerEvent;
0187 
0188   HLTConfigProvider hltConfig_;
0189   edm::Handle<edm::TriggerResults> HLTR;
0190   std::string getTriggerName(std::string partialName);
0191 };
0192 
0193 BPHMonitor::BPHMonitor(const edm::ParameterSet& iConfig)
0194     : folderName_(iConfig.getParameter<std::string>("FolderName")),
0195       requireValidHLTPaths_(iConfig.getParameter<bool>("requireValidHLTPaths")),
0196       hltPathsAreValid_(false),
0197       muoInputTag_(iConfig.getParameter<edm::InputTag>("muons")),
0198       bsInputTag_(iConfig.getParameter<edm::InputTag>("beamSpot")),
0199       trInputTag_(iConfig.getParameter<edm::InputTag>("tracks")),
0200       phInputTag_(iConfig.getParameter<edm::InputTag>("photons")),
0201       vtxInputTag_(iConfig.getParameter<edm::InputTag>("offlinePVs")),
0202       muoToken_(mayConsume<reco::MuonCollection>(muoInputTag_)),
0203       bsToken_(mayConsume<reco::BeamSpot>(bsInputTag_)),
0204       trToken_(mayConsume<reco::TrackCollection>(trInputTag_)),
0205       phToken_(mayConsume<reco::PhotonCollection>(phInputTag_)),
0206       vtxToken_(mayConsume<reco::VertexCollection>(vtxInputTag_)),
0207       pt_variable_binning_(
0208           iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double>>("ptBinning")),
0209       dMu_pt_variable_binning_(
0210           iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double>>("dMuPtBinning")),
0211       prob_variable_binning_(
0212           iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double>>("probBinning")),
0213       phi_binning_(getHistoPSet(
0214           iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("phiPSet"))),
0215       eta_binning_(getHistoPSet(
0216           iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("etaPSet"))),
0217       d0_binning_(
0218           getHistoPSet(iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("d0PSet"))),
0219       z0_binning_(
0220           getHistoPSet(iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("z0PSet"))),
0221       dR_binning_(
0222           getHistoPSet(iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("dRPSet"))),
0223       mass_binning_(getHistoPSet(
0224           iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("massPSet"))),
0225       Bmass_binning_(getHistoPSet(
0226           iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("BmassPSet"))),
0227       dca_binning_(getHistoPSet(
0228           iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("dcaPSet"))),
0229       ds_binning_(
0230           getHistoPSet(iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("dsPSet"))),
0231       cos_binning_(getHistoPSet(
0232           iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<edm::ParameterSet>("cosPSet"))),
0233       num_genTriggerEventFlag_(new GenericTriggerEventFlag(
0234           iConfig.getParameter<edm::ParameterSet>("numGenericTriggerEventPSet"), consumesCollector(), *this)),
0235       den_genTriggerEventFlag_(new GenericTriggerEventFlag(
0236           iConfig.getParameter<edm::ParameterSet>("denGenericTriggerEventPSet"), consumesCollector(), *this)),
0237       hltPrescale_(new HLTPrescaleProvider(iConfig, consumesCollector(), *this)),
0238       muoSelection_(iConfig.getParameter<std::string>("muoSelection")),
0239       muoSelection_ref(iConfig.getParameter<std::string>("muoSelection_ref")),
0240       muoSelection_tag(iConfig.getParameter<std::string>("muoSelection_tag")),
0241       muoSelection_probe(iConfig.getParameter<std::string>("muoSelection_probe")),
0242       nmuons_(iConfig.getParameter<int>("nmuons")),
0243       tnp_(iConfig.getParameter<bool>("tnp")),
0244       L3_(iConfig.getParameter<int>("L3")),
0245       ptCut_(iConfig.getParameter<int>("ptCut")),
0246       displaced_(iConfig.getParameter<int>("displaced")),
0247       trOrMu_(iConfig.getParameter<int>("trOrMu")),
0248       Jpsi_(iConfig.getParameter<int>("Jpsi")),
0249       Upsilon_(iConfig.getParameter<int>("Upsilon"))  // if ==1 path with Upsilon constraint
0250       ,
0251       enum_(iConfig.getParameter<int>("enum")),
0252       seagull_(iConfig.getParameter<int>("seagull")),
0253       maxmass_(iConfig.getParameter<double>("maxmass")),
0254       minmass_(iConfig.getParameter<double>("minmass")),
0255       maxmassJpsi(iConfig.getParameter<double>("maxmassJpsi")),
0256       minmassJpsi(iConfig.getParameter<double>("minmassJpsi")),
0257       maxmassUpsilon(iConfig.getParameter<double>("maxmassUpsilon")),
0258       minmassUpsilon(iConfig.getParameter<double>("minmassUpsilon")),
0259       maxmassTkTk(iConfig.getParameter<double>("maxmassTkTk")),
0260       minmassTkTk(iConfig.getParameter<double>("minmassTkTk")),
0261       maxmassJpsiTk(iConfig.getParameter<double>("maxmassJpsiTk")),
0262       minmassJpsiTk(iConfig.getParameter<double>("minmassJpsiTk")),
0263       kaon_mass(iConfig.getParameter<double>("kaon_mass")),
0264       mu_mass(iConfig.getParameter<double>("mu_mass")),
0265       min_dR(iConfig.getParameter<double>("min_dR")),
0266       max_dR(iConfig.getParameter<double>("max_dR")),
0267       minprob(iConfig.getParameter<double>("minprob")),
0268       mincos(iConfig.getParameter<double>("mincos")),
0269       minDS(iConfig.getParameter<double>("minDS")),
0270       hltInputTag_1(iConfig.getParameter<edm::InputTag>("hltTriggerSummaryAOD")),
0271       hltInputTag_(mayConsume<trigger::TriggerEvent>(iConfig.getParameter<edm::InputTag>("hltTriggerSummaryAOD"))),
0272       hltpaths_num(iConfig.getParameter<edm::ParameterSet>("numGenericTriggerEventPSet")
0273                        .getParameter<std::vector<std::string>>("hltPaths")),
0274       hltpaths_den(iConfig.getParameter<edm::ParameterSet>("denGenericTriggerEventPSet")
0275                        .getParameter<std::vector<std::string>>("hltPaths")),
0276       trSelection_(iConfig.getParameter<std::string>("muoSelection")),
0277       trSelection_ref(iConfig.getParameter<std::string>("trSelection_ref")),
0278       DMSelection_ref(iConfig.getParameter<std::string>("DMSelection_ref")) {
0279   if (!tnp_) {
0280     magneticFieldToken_ = esConsumes<MagneticField, IdealMagneticFieldRecord>();
0281   }
0282 }
0283 
0284 BPHMonitor::~BPHMonitor() throw() {
0285   if (num_genTriggerEventFlag_) {
0286     num_genTriggerEventFlag_.reset();
0287   }
0288   if (den_genTriggerEventFlag_) {
0289     den_genTriggerEventFlag_.reset();
0290   }
0291 
0292   delete hltPrescale_;
0293 }
0294 
0295 void BPHMonitor::bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& iRun, edm::EventSetup const& iSetup) {
0296   // Initialize the GenericTriggerEventFlag
0297   if (num_genTriggerEventFlag_ && num_genTriggerEventFlag_->on()) {
0298     num_genTriggerEventFlag_->initRun(iRun, iSetup);
0299   }
0300   if (den_genTriggerEventFlag_ && den_genTriggerEventFlag_->on()) {
0301     den_genTriggerEventFlag_->initRun(iRun, iSetup);
0302   }
0303 
0304   // check if every HLT path specified in numerator and denominator has a valid match in the HLT Menu
0305   hltPathsAreValid_ = (num_genTriggerEventFlag_ && den_genTriggerEventFlag_ && num_genTriggerEventFlag_->on() &&
0306                        den_genTriggerEventFlag_->on() && num_genTriggerEventFlag_->allHLTPathsAreValid() &&
0307                        den_genTriggerEventFlag_->allHLTPathsAreValid());
0308 
0309   // if valid HLT paths are required,
0310   // create DQM outputs only if all paths are valid
0311   if (requireValidHLTPaths_ and (not hltPathsAreValid_)) {
0312     return;
0313   }
0314 
0315   std::string histname, histtitle, istnp, trMuPh;
0316 
0317   bool Ph_ = false;
0318   if (enum_ == 7)
0319     Ph_ = true;
0320   if (tnp_)
0321     istnp = "Tag_and_Probe/";
0322   else
0323     istnp = "";
0324 
0325   std::string currentFolder = folderName_ + istnp;
0326   ibooker.setCurrentFolder(currentFolder);
0327 
0328   if (trOrMu_)
0329     trMuPh = "tr";
0330   else if (Ph_)
0331     trMuPh = "ph";
0332   else
0333     trMuPh = "mu";
0334 
0335   if (enum_ == 7 || enum_ == 1 || enum_ == 9 || enum_ == 10) {
0336     histname = trMuPh + "Pt";
0337     histtitle = trMuPh + "_P_{t}";
0338     bookME(ibooker, muPt_, histname, histtitle, pt_variable_binning_);
0339     setMETitle(muPt_, trMuPh + "_Pt[GeV]", "events / 1 GeV");
0340 
0341     histname = trMuPh + "Phi";
0342     histtitle = trMuPh + "Phi";
0343     bookME(ibooker, muPhi_, histname, histtitle, phi_binning_.nbins, phi_binning_.xmin, phi_binning_.xmax);
0344     setMETitle(muPhi_, trMuPh + "_#phi", "events / 0.1 rad");
0345 
0346     histname = trMuPh + "Eta";
0347     histtitle = trMuPh + "_Eta";
0348     bookME(ibooker, muEta_, histname, histtitle, eta_binning_.nbins, eta_binning_.xmin, eta_binning_.xmax);
0349     setMETitle(muEta_, trMuPh + "_#eta", "events / 0.2");
0350 
0351     if (enum_ == 9) {
0352       histname = "BMass";
0353       histtitle = "BMass";
0354       bookME(ibooker, BMass_, histname, histtitle, Bmass_binning_.nbins, mass_binning_.xmin, mass_binning_.xmax);
0355       setMETitle(BMass_, "B_#mass", "events /");
0356     }
0357   } else {
0358     if (enum_ != 8) {
0359       histname = trMuPh + "1Pt";
0360       histtitle = trMuPh + "1_P_{t}";
0361       bookME(ibooker, mu1Pt_, histname, histtitle, pt_variable_binning_);
0362       setMETitle(mu1Pt_, trMuPh + "_Pt[GeV]", "events / 1 GeV");
0363 
0364       histname = trMuPh + "1Phi";
0365       histtitle = trMuPh + "1Phi";
0366       bookME(ibooker, mu1Phi_, histname, histtitle, phi_binning_.nbins, phi_binning_.xmin, phi_binning_.xmax);
0367       setMETitle(mu1Phi_, trMuPh + "_#phi", "events / 0.1 rad");
0368 
0369       histname = trMuPh + "1Eta";
0370       histtitle = trMuPh + "1_Eta";
0371       bookME(ibooker, mu1Eta_, histname, histtitle, eta_binning_.nbins, eta_binning_.xmin, eta_binning_.xmax);
0372       setMETitle(mu1Eta_, trMuPh + "_#eta", "events / 0.2");
0373 
0374       histname = trMuPh + "2Pt";
0375       histtitle = trMuPh + "2_P_{t}";
0376       bookME(ibooker, mu2Pt_, histname, histtitle, pt_variable_binning_);
0377       setMETitle(mu2Pt_, trMuPh + "_Pt[GeV]", "events / 1 GeV");
0378 
0379       histname = trMuPh + "2Phi";
0380       histtitle = trMuPh + "2Phi";
0381       bookME(ibooker, mu2Phi_, histname, histtitle, phi_binning_.nbins, phi_binning_.xmin, phi_binning_.xmax);
0382       setMETitle(mu2Phi_, trMuPh + "_#phi", "events / 0.1 rad");
0383 
0384       histname = trMuPh + "2Eta";
0385       histtitle = trMuPh + "2_Eta";
0386       bookME(ibooker, mu2Eta_, histname, histtitle, eta_binning_.nbins, eta_binning_.xmin, eta_binning_.xmax);
0387       setMETitle(mu2Eta_, trMuPh + "_#eta", "events / 0.2");
0388       if (enum_ == 11) {
0389         histname = "BMass";
0390         histtitle = "BMass";
0391         bookME(ibooker, BMass_, histname, histtitle, Bmass_binning_.nbins, mass_binning_.xmin, mass_binning_.xmax);
0392         setMETitle(BMass_, "B_#mass", "events /");
0393       }
0394     }
0395     if (enum_ == 6) {
0396       histname = trMuPh + "3Eta";
0397       histtitle = trMuPh + "3Eta";
0398       bookME(ibooker, mu3Eta_, histname, histtitle, eta_binning_.nbins, eta_binning_.xmin, eta_binning_.xmax);
0399       setMETitle(mu3Eta_, trMuPh + "3#eta", "events / 0.2");
0400 
0401       histname = trMuPh + "3Pt";
0402       histtitle = trMuPh + "3_P_{t}";
0403       bookME(ibooker, mu3Pt_, histname, histtitle, pt_variable_binning_);
0404       setMETitle(mu3Pt_, trMuPh + "3_Pt[GeV]", "events / 1 GeV");
0405 
0406       histname = trMuPh + "3Phi";
0407       histtitle = trMuPh + "3Phi";
0408       bookME(ibooker, mu3Phi_, histname, histtitle, phi_binning_.nbins, phi_binning_.xmin, phi_binning_.xmax);
0409       setMETitle(mu3Phi_, trMuPh + "3_#phi", "events / 0.1 rad");
0410 
0411     } else if (enum_ == 2 || enum_ == 4 || enum_ == 5 || enum_ == 8) {
0412       histname = "DiMuEta";
0413       histtitle = "DiMuEta";
0414       bookME(ibooker, DiMuEta_, histname, histtitle, eta_binning_.nbins, eta_binning_.xmin, eta_binning_.xmax);
0415       setMETitle(DiMuEta_, "DiMu#eta", "events / 0.2");
0416 
0417       histname = "DiMuPt";
0418       histtitle = "DiMu_P_{t}";
0419       bookME(ibooker, DiMuPt_, histname, histtitle, dMu_pt_variable_binning_);
0420       setMETitle(DiMuPt_, "DiMu_Pt[GeV]", "events / 1 GeV");
0421 
0422       histname = "DiMuPhi";
0423       histtitle = "DiMuPhi";
0424       bookME(ibooker, DiMuPhi_, histname, histtitle, phi_binning_.nbins, phi_binning_.xmin, phi_binning_.xmax);
0425       setMETitle(DiMuPhi_, "DiMu_#phi", "events / 0.1 rad");
0426 
0427       if (enum_ == 4 || enum_ == 5) {
0428         histname = "DiMudR";
0429         histtitle = "DiMudR";
0430         bookME(ibooker, DiMudR_, histname, histtitle, dR_binning_.nbins, dR_binning_.xmin, dR_binning_.xmax);
0431         setMETitle(DiMudR_, "DiMu_#dR", "events /");
0432 
0433         if (enum_ == 4) {
0434           histname = "DiMuMass";
0435           histtitle = "DiMuMass";
0436           bookME(ibooker, DiMuMass_, histname, histtitle, mass_binning_.nbins, mass_binning_.xmin, mass_binning_.xmax);
0437           setMETitle(DiMuMass_, "DiMu_#mass", "events /");
0438         }
0439       } else if (enum_ == 8) {
0440         histname = "DiMuProb";
0441         histtitle = "DiMuProb";
0442         bookME(ibooker, DiMuProb_, histname, histtitle, prob_variable_binning_);
0443         setMETitle(DiMuProb_, "DiMu_#prob", "events /");
0444 
0445         histname = "DiMuPVcos";
0446         histtitle = "DiMuPVcos";
0447         bookME(ibooker, DiMuPVcos_, histname, histtitle, cos_binning_.nbins, cos_binning_.xmin, cos_binning_.xmax);
0448         setMETitle(DiMuPVcos_, "DiMu_#cosPV", "events /");
0449 
0450         histname = "DiMuDS";
0451         histtitle = "DiMuDS";
0452         bookME(ibooker, DiMuDS_, histname, histtitle, ds_binning_.nbins, ds_binning_.xmin, ds_binning_.xmax);
0453         setMETitle(DiMuDS_, "DiMu_#ds", "events /");
0454 
0455         histname = "DiMuDCA";
0456         histtitle = "DiMuDCA";
0457         bookME(ibooker, DiMuDCA_, histname, histtitle, dca_binning_.nbins, dca_binning_.xmin, dca_binning_.xmax);
0458         setMETitle(DiMuDCA_, "DiMu_#dca", "events /");
0459       }
0460     }  // if (enum_ == 2 || enum_ == 4 || enum_ == 5 || enum_ == 8)
0461   }
0462 
0463   bool changed = true;
0464 
0465   hltPrescale_->init(iRun, iSetup, "HLT", changed);
0466   hltConfig_ = hltPrescale_->hltConfigProvider();
0467 }
0468 
0469 void BPHMonitor::analyze(edm::Event const& iEvent, edm::EventSetup const& iSetup) {
0470   // if valid HLT paths are required,
0471   // analyze event only if all paths are valid
0472   if (requireValidHLTPaths_ and (not hltPathsAreValid_)) {
0473     return;
0474   }
0475 
0476   edm::Handle<reco::BeamSpot> beamSpot;
0477   iEvent.getByToken(bsToken_, beamSpot);
0478   if (!beamSpot.isValid()) {
0479     return;
0480   }
0481 
0482   edm::Handle<reco::MuonCollection> muoHandle;
0483   iEvent.getByToken(muoToken_, muoHandle);
0484   if (!muoHandle.isValid()) {
0485     return;
0486   }
0487 
0488   edm::Handle<reco::TrackCollection> trHandle;
0489   iEvent.getByToken(trToken_, trHandle);
0490   if (!trHandle.isValid()) {
0491     return;
0492   }
0493 
0494   edm::Handle<reco::PhotonCollection> phHandle;
0495   iEvent.getByToken(phToken_, phHandle);
0496 
0497   edm::Handle<edm::TriggerResults> handleTriggerTrigRes;
0498 
0499   const std::string& hltpath = getTriggerName(hltpaths_den[0]);
0500   const std::string& hltpath1 = getTriggerName(hltpaths_num[0]);
0501 
0502   double PrescaleWeight = 1.0;
0503   if (den_genTriggerEventFlag_->on() && den_genTriggerEventFlag_->accept(iEvent, iSetup) &&
0504       num_genTriggerEventFlag_->on() && num_genTriggerEventFlag_->accept(iEvent, iSetup))
0505     PrescaleWeight = Prescale(hltpath1, hltpath, iEvent, iSetup, hltPrescale_);
0506 
0507   if (tnp_ > 0) {  //TnP method
0508 
0509     if (den_genTriggerEventFlag_->on() && !den_genTriggerEventFlag_->accept(iEvent, iSetup))
0510       return;
0511     iEvent.getByToken(hltInputTag_, handleTriggerEvent);
0512     if (handleTriggerEvent->sizeFilters() == 0)
0513       return;
0514 
0515     std::vector<reco::Muon> tagMuons;
0516     for (auto const& m : *muoHandle) {  // applying tag selection
0517       if (!matchToTrigger(hltpath, m))
0518         continue;
0519       if (muoSelection_ref(m))
0520         tagMuons.push_back(m);
0521     }
0522 
0523     for (int i = 0; i < int(tagMuons.size()); i++) {
0524       for (auto const& m : *muoHandle) {
0525         if (!matchToTrigger(hltpath, m))
0526           continue;
0527         if ((tagMuons[i].pt() == m.pt()))
0528           continue;  //not the same
0529         if ((tagMuons[i].p4() + m.p4()).M() > minmass_ &&
0530             (tagMuons[i].p4() + m.p4()).M() < maxmass_) {  //near to J/psi mass
0531           muPhi_.denominator->Fill(m.phi());
0532           muEta_.denominator->Fill(m.eta());
0533           muPt_.denominator->Fill(m.pt());
0534           if (muoSelection_(m) && num_genTriggerEventFlag_->on() && num_genTriggerEventFlag_->accept(iEvent, iSetup)) {
0535             muPhi_.numerator->Fill(m.phi(), PrescaleWeight);
0536             muEta_.numerator->Fill(m.eta(), PrescaleWeight);
0537             muPt_.numerator->Fill(m.pt(), PrescaleWeight);
0538           }
0539         }
0540       }
0541     }
0542 
0543   } else {  // reference method
0544 
0545     if (den_genTriggerEventFlag_->on() && (!den_genTriggerEventFlag_->accept(iEvent, iSetup)))
0546       return;
0547 
0548     iEvent.getByToken(hltInputTag_, handleTriggerEvent);
0549     if (handleTriggerEvent->sizeFilters() == 0)
0550       return;
0551 
0552     for (auto const& m : *muoHandle) {
0553       if (!muoSelection_ref(m))
0554         continue;
0555       if (!matchToTrigger(hltpath, m))
0556         continue;
0557 
0558       for (auto const& m1 : *muoHandle) {
0559         if (!(m1.pt() > m.pt()))
0560           continue;
0561         if (ptCut_) {
0562           if (!muoSelection_(m1))
0563             continue;
0564         } else if (!muoSelection_ref(m1))
0565           continue;
0566         if (!matchToTrigger(hltpath, m1))
0567           continue;
0568 
0569         if (enum_ != 10) {
0570           if (!DMSelection_ref(m1.p4() + m.p4()))
0571             continue;
0572           if (m.charge() * m1.charge() > 0)
0573             continue;
0574         }
0575 
0576         // dimuon vertex reconstruction
0577         MagneticField const& magneticField = iSetup.getData(magneticFieldToken_);
0578         const reco::BeamSpot& vertexBeamSpot = *beamSpot;
0579         std::vector<reco::TransientTrack> j_tks;
0580         reco::TransientTrack mu1TT(m.track(), &magneticField);
0581         reco::TransientTrack mu2TT(m1.track(), &magneticField);
0582         j_tks.push_back(mu1TT);
0583         j_tks.push_back(mu2TT);
0584         KalmanVertexFitter jkvf;
0585         TransientVertex jtv = jkvf.vertex(j_tks);
0586         if (!jtv.isValid())
0587           continue;
0588         reco::Vertex jpsivertex = jtv;
0589         float dimuonCL = 0;
0590         if ((jpsivertex.chi2() >= 0) && (jpsivertex.ndof() > 0))
0591           dimuonCL = TMath::Prob(jpsivertex.chi2(), jpsivertex.ndof());
0592         math::XYZVector jpperp(m.px() + m1.px(), m.py() + m1.py(), 0.);
0593         GlobalPoint jVertex = jtv.position();
0594         GlobalError jerr = jtv.positionError();
0595         GlobalPoint displacementFromBeamspotJpsi(
0596             -1 * ((vertexBeamSpot.x0() - jVertex.x()) + (jVertex.z() - vertexBeamSpot.z0()) * vertexBeamSpot.dxdz()),
0597             -1 * ((vertexBeamSpot.y0() - jVertex.y()) + (jVertex.z() - vertexBeamSpot.z0()) * vertexBeamSpot.dydz()),
0598             0);
0599         reco::Vertex::Point vperpj(displacementFromBeamspotJpsi.x(), displacementFromBeamspotJpsi.y(), 0.);
0600         float jpsi_cos = vperpj.Dot(jpperp) / (vperpj.R() * jpperp.R());
0601         TrajectoryStateClosestToPoint mu1TS = mu1TT.impactPointTSCP();
0602         TrajectoryStateClosestToPoint mu2TS = mu2TT.impactPointTSCP();
0603         ClosestApproachInRPhi cApp;
0604         if (mu1TS.isValid() && mu2TS.isValid()) {
0605           if (!cApp.calculate(mu1TS.theState(), mu2TS.theState()))
0606             continue;
0607         } else
0608           continue;
0609         double DiMuMass = (m1.p4() + m.p4()).M();
0610 
0611         switch (
0612             enum_) {  // enum_ = 1...9, represents different sets of variables for different paths, we want to have different hists for different paths
0613 
0614           case 1:
0615 
0616             tnp_ = true;  // already filled hists for tnp method
0617             [[fallthrough]];
0618           case 2:
0619 
0620             if ((Jpsi_) && (!Upsilon_))
0621               if (DiMuMass > maxmassJpsi || DiMuMass < minmassJpsi)
0622                 continue;
0623             if ((!Jpsi_) && (Upsilon_))
0624               if (DiMuMass > maxmassUpsilon || DiMuMass < minmassUpsilon)
0625                 continue;
0626             if (dimuonCL < minprob)
0627               continue;
0628 
0629             mu1Phi_.denominator->Fill(m.phi());
0630             mu1Eta_.denominator->Fill(m.eta());
0631             mu1Pt_.denominator->Fill(m.pt());
0632             mu2Phi_.denominator->Fill(m1.phi());
0633             mu2Eta_.denominator->Fill(m1.eta());
0634             mu2Pt_.denominator->Fill(m1.pt());
0635             DiMuPt_.denominator->Fill((m1.p4() + m.p4()).Pt());
0636             DiMuEta_.denominator->Fill((m1.p4() + m.p4()).Eta());
0637             DiMuPhi_.denominator->Fill((m1.p4() + m.p4()).Phi());
0638 
0639             if (num_genTriggerEventFlag_->on() && num_genTriggerEventFlag_->accept(iEvent, iSetup)) {
0640               if (!matchToTrigger(hltpath1, m1))
0641                 continue;
0642               if (!matchToTrigger(hltpath1, m))
0643                 continue;
0644               mu1Phi_.numerator->Fill(m.phi(), PrescaleWeight);
0645               mu1Eta_.numerator->Fill(m.eta(), PrescaleWeight);
0646               mu1Pt_.numerator->Fill(m.pt(), PrescaleWeight);
0647               mu2Phi_.numerator->Fill(m1.phi(), PrescaleWeight);
0648               mu2Eta_.numerator->Fill(m1.eta(), PrescaleWeight);
0649               mu2Pt_.numerator->Fill(m1.pt(), PrescaleWeight);
0650               DiMuPt_.numerator->Fill((m1.p4() + m.p4()).Pt(), PrescaleWeight);
0651               DiMuEta_.numerator->Fill((m1.p4() + m.p4()).Eta(), PrescaleWeight);
0652               DiMuPhi_.numerator->Fill((m1.p4() + m.p4()).Phi(), PrescaleWeight);
0653             }
0654 
0655             break;
0656 
0657           case 3:
0658 
0659             if ((Jpsi_) && (!Upsilon_))
0660               if (DiMuMass > maxmassJpsi || DiMuMass < minmassJpsi)
0661                 continue;
0662             if ((!Jpsi_) && (Upsilon_))
0663               if (DiMuMass > maxmassUpsilon || DiMuMass < minmassUpsilon)
0664                 continue;
0665             if (dimuonCL < minprob)
0666               continue;
0667 
0668             mu1Phi_.denominator->Fill(m.phi());
0669             mu1Eta_.denominator->Fill(m.eta());
0670             mu1Pt_.denominator->Fill(m.pt());
0671             mu2Phi_.denominator->Fill(m1.phi());
0672             mu2Eta_.denominator->Fill(m1.eta());
0673             mu2Pt_.denominator->Fill(m1.pt());
0674 
0675             if (num_genTriggerEventFlag_->on() && num_genTriggerEventFlag_->accept(iEvent, iSetup)) {
0676               if (!matchToTrigger(hltpath1, m1))
0677                 continue;
0678               if (!matchToTrigger(hltpath1, m))
0679                 continue;
0680               mu1Phi_.numerator->Fill(m.phi(), PrescaleWeight);
0681               mu1Eta_.numerator->Fill(m.eta(), PrescaleWeight);
0682               mu1Pt_.numerator->Fill(m.pt(), PrescaleWeight);
0683               mu2Phi_.numerator->Fill(m1.phi(), PrescaleWeight);
0684               mu2Eta_.numerator->Fill(m1.eta(), PrescaleWeight);
0685               mu2Pt_.numerator->Fill(m1.pt(), PrescaleWeight);
0686             }
0687 
0688             break;
0689 
0690           case 4:
0691 
0692             if (dimuonCL < minprob)
0693               continue;
0694 
0695             // fill mass plots without selecting mass region
0696             DiMuMass_.denominator->Fill(DiMuMass);
0697             if (num_genTriggerEventFlag_->on() && num_genTriggerEventFlag_->accept(iEvent, iSetup) &&
0698                 !(seagull_ && m.charge() * deltaPhi(m.phi(), m1.phi()) > 0) && matchToTrigger(hltpath1, m1) &&
0699                 matchToTrigger(hltpath1, m))
0700               DiMuMass_.numerator->Fill(DiMuMass, PrescaleWeight);
0701 
0702             if ((Jpsi_) && (!Upsilon_))
0703               if (DiMuMass > maxmassJpsi || DiMuMass < minmassJpsi)
0704                 continue;
0705             if ((!Jpsi_) && (Upsilon_))
0706               if (DiMuMass > maxmassUpsilon || DiMuMass < minmassUpsilon)
0707                 continue;
0708 
0709             mu1Phi_.denominator->Fill(m.phi());
0710             mu1Eta_.denominator->Fill(m.eta());
0711             mu1Pt_.denominator->Fill(m.pt());
0712             mu2Phi_.denominator->Fill(m1.phi());
0713             mu2Eta_.denominator->Fill(m1.eta());
0714             mu2Pt_.denominator->Fill(m1.pt());
0715             DiMuPt_.denominator->Fill((m1.p4() + m.p4()).Pt());
0716             DiMuEta_.denominator->Fill((m1.p4() + m.p4()).Eta());
0717             DiMuPhi_.denominator->Fill((m1.p4() + m.p4()).Phi());
0718             DiMudR_.denominator->Fill(reco::deltaR(m, m1));
0719 
0720             if (num_genTriggerEventFlag_->on() && num_genTriggerEventFlag_->accept(iEvent, iSetup)) {
0721               if (seagull_ && m.charge() * deltaPhi(m.phi(), m1.phi()) > 0)
0722                 continue;
0723               if (!matchToTrigger(hltpath1, m1))
0724                 continue;
0725               if (!matchToTrigger(hltpath1, m))
0726                 continue;
0727               mu1Phi_.numerator->Fill(m.phi(), PrescaleWeight);
0728               mu1Eta_.numerator->Fill(m.eta(), PrescaleWeight);
0729               mu1Pt_.numerator->Fill(m.pt(), PrescaleWeight);
0730               mu2Phi_.numerator->Fill(m1.phi(), PrescaleWeight);
0731               mu2Eta_.numerator->Fill(m1.eta(), PrescaleWeight);
0732               mu2Pt_.numerator->Fill(m1.pt(), PrescaleWeight);
0733               DiMuPt_.numerator->Fill((m1.p4() + m.p4()).Pt(), PrescaleWeight);
0734               DiMuEta_.numerator->Fill((m1.p4() + m.p4()).Eta(), PrescaleWeight);
0735               DiMuPhi_.numerator->Fill((m1.p4() + m.p4()).Phi(), PrescaleWeight);
0736               DiMudR_.numerator->Fill(reco::deltaR(m, m1), PrescaleWeight);
0737             }
0738 
0739             break;
0740 
0741           case 5:
0742 
0743             if (dimuonCL < minprob)
0744               continue;
0745             if ((Jpsi_) && (!Upsilon_))
0746               if (DiMuMass > maxmassJpsi || DiMuMass < minmassJpsi)
0747                 continue;
0748             if ((!Jpsi_) && (Upsilon_))
0749               if (DiMuMass > maxmassUpsilon || DiMuMass < minmassUpsilon)
0750                 continue;
0751 
0752             mu1Phi_.denominator->Fill(m.phi());
0753             mu1Eta_.denominator->Fill(m.eta());
0754             mu1Pt_.denominator->Fill(m.pt());
0755             mu2Phi_.denominator->Fill(m1.phi());
0756             mu2Eta_.denominator->Fill(m1.eta());
0757             mu2Pt_.denominator->Fill(m1.pt());
0758             DiMuPt_.denominator->Fill((m1.p4() + m.p4()).Pt());
0759             DiMuEta_.denominator->Fill((m1.p4() + m.p4()).Eta());
0760             DiMuPhi_.denominator->Fill((m1.p4() + m.p4()).Phi());
0761             DiMudR_.denominator->Fill(reco::deltaR(m, m1));
0762 
0763             if (num_genTriggerEventFlag_->on() && num_genTriggerEventFlag_->accept(iEvent, iSetup)) {
0764               if (seagull_ && m.charge() * deltaPhi(m.phi(), m1.phi()) > 0)
0765                 continue;
0766               if (!matchToTrigger(hltpath1, m1))
0767                 continue;
0768               if (!matchToTrigger(hltpath1, m))
0769                 continue;
0770               mu1Phi_.numerator->Fill(m.phi(), PrescaleWeight);
0771               mu1Eta_.numerator->Fill(m.eta(), PrescaleWeight);
0772               mu1Pt_.numerator->Fill(m.pt(), PrescaleWeight);
0773               mu2Phi_.numerator->Fill(m1.phi(), PrescaleWeight);
0774               mu2Eta_.numerator->Fill(m1.eta(), PrescaleWeight);
0775               mu2Pt_.numerator->Fill(m1.pt(), PrescaleWeight);
0776               DiMuPt_.numerator->Fill((m1.p4() + m.p4()).Pt(), PrescaleWeight);
0777               DiMuEta_.numerator->Fill((m1.p4() + m.p4()).Eta(), PrescaleWeight);
0778               DiMuPhi_.numerator->Fill((m1.p4() + m.p4()).Phi(), PrescaleWeight);
0779               DiMudR_.numerator->Fill(reco::deltaR(m, m1), PrescaleWeight);
0780             }
0781 
0782             break;
0783 
0784           case 6:
0785 
0786             if (dimuonCL < minprob)
0787               continue;
0788             if ((Jpsi_) && (!Upsilon_))
0789               if (DiMuMass > maxmassJpsi || DiMuMass < minmassJpsi)
0790                 continue;
0791             if ((!Jpsi_) && (Upsilon_))
0792               if (DiMuMass > maxmassUpsilon || DiMuMass < minmassUpsilon)
0793                 continue;
0794 
0795             for (auto const& m2 : *muoHandle) {
0796               if (m2.pt() == m.pt())
0797                 continue;  // remove duplicates but do not introduce ordering
0798               if (m2.pt() == m1.pt())
0799                 continue;  // -> m2 will be the non-resonant and non-vertexing muon in the triplet
0800               if (!matchToTrigger(hltpath, m2))
0801                 continue;
0802 
0803               mu1Phi_.denominator->Fill(m.phi());
0804               mu1Eta_.denominator->Fill(m.eta());
0805               mu1Pt_.denominator->Fill(m.pt());
0806               mu2Phi_.denominator->Fill(m1.phi());
0807               mu2Eta_.denominator->Fill(m1.eta());
0808               mu2Pt_.denominator->Fill(m1.pt());
0809               mu3Phi_.denominator->Fill(m2.phi());
0810               mu3Eta_.denominator->Fill(m2.eta());
0811               mu3Pt_.denominator->Fill(m2.pt());
0812 
0813               if (num_genTriggerEventFlag_->on() && num_genTriggerEventFlag_->accept(iEvent, iSetup)) {
0814                 if (!matchToTrigger(hltpath1, m1))
0815                   continue;
0816                 if (!matchToTrigger(hltpath1, m))
0817                   continue;
0818                 if (!matchToTrigger(hltpath1, m2))
0819                   continue;
0820                 mu1Phi_.numerator->Fill(m.phi(), PrescaleWeight);
0821                 mu1Eta_.numerator->Fill(m.eta(), PrescaleWeight);
0822                 mu1Pt_.numerator->Fill(m.pt(), PrescaleWeight);
0823                 mu2Phi_.numerator->Fill(m1.phi(), PrescaleWeight);
0824                 mu2Eta_.numerator->Fill(m1.eta(), PrescaleWeight);
0825                 mu2Pt_.numerator->Fill(m1.pt(), PrescaleWeight);
0826                 mu3Phi_.numerator->Fill(m2.phi(), PrescaleWeight);
0827                 mu3Eta_.numerator->Fill(m2.eta(), PrescaleWeight);
0828                 mu3Pt_.numerator->Fill(m2.pt(), PrescaleWeight);
0829               }
0830             }
0831 
0832             break;
0833 
0834           case 7:
0835 
0836             if (phHandle.isValid()) {
0837               if (!phHandle->empty())
0838                 for (auto const& p : *phHandle) {
0839                   if (!matchToTrigger(hltpath, p))
0840                     continue;
0841 
0842                   phPhi_.denominator->Fill(p.phi());
0843                   phEta_.denominator->Fill(p.eta());
0844                   phPt_.denominator->Fill(p.pt());
0845 
0846                   if (num_genTriggerEventFlag_->on() && num_genTriggerEventFlag_->accept(iEvent, iSetup)) {
0847                     if (!matchToTrigger(hltpath1, p))
0848                       continue;
0849                     if (!matchToTrigger(hltpath1, m))
0850                       continue;
0851                     if (!matchToTrigger(hltpath1, m1))
0852                       continue;
0853                     phPhi_.numerator->Fill(p.phi(), PrescaleWeight);
0854                     phEta_.numerator->Fill(p.eta(), PrescaleWeight);
0855                     phPt_.numerator->Fill(p.pt(), PrescaleWeight);
0856                   }
0857                 }
0858             } else {
0859               // if Handle is not valid, because the InputTag has been mis-configured, then skip the event
0860               if (!phInputTag_.label().empty())
0861                 return;
0862             }
0863 
0864             break;
0865 
0866           case 8:  //vtx monitoring, filling probability, DS, DCA, cos of pointing angle to the PV, eta, pT of dimuon
0867 
0868             if (displaced_)
0869               if ((displacementFromBeamspotJpsi.perp() / sqrt(jerr.rerr(displacementFromBeamspotJpsi))) < minDS)
0870                 continue;
0871             if ((Jpsi_) && (!Upsilon_))
0872               if (DiMuMass > maxmassJpsi || DiMuMass < minmassJpsi)
0873                 continue;
0874             if ((!Jpsi_) && (Upsilon_))
0875               if (DiMuMass > maxmassUpsilon || DiMuMass < minmassUpsilon)
0876                 continue;
0877 
0878             // fill vtx-prob plots before selecting on this variable
0879             DiMuProb_.denominator->Fill(dimuonCL);
0880             if (num_genTriggerEventFlag_->on() && num_genTriggerEventFlag_->accept(iEvent, iSetup) &&
0881                 matchToTrigger(hltpath1, m1) && matchToTrigger(hltpath1, m))
0882               DiMuProb_.numerator->Fill(dimuonCL, PrescaleWeight);
0883 
0884             if (dimuonCL < minprob)
0885               continue;
0886 
0887             DiMuDS_.denominator->Fill(displacementFromBeamspotJpsi.perp() /
0888                                       sqrt(jerr.rerr(displacementFromBeamspotJpsi)));
0889             DiMuPVcos_.denominator->Fill(jpsi_cos);
0890             DiMuPt_.denominator->Fill((m1.p4() + m.p4()).Pt());
0891             DiMuEta_.denominator->Fill((m1.p4() + m.p4()).Eta());
0892             DiMuPhi_.denominator->Fill((m1.p4() + m.p4()).Phi());
0893             DiMuDCA_.denominator->Fill(cApp.distance());
0894 
0895             if (num_genTriggerEventFlag_->on() && num_genTriggerEventFlag_->accept(iEvent, iSetup)) {
0896               if (!matchToTrigger(hltpath1, m1))
0897                 continue;
0898               if (!matchToTrigger(hltpath1, m))
0899                 continue;
0900               DiMuDS_.numerator->Fill(
0901                   displacementFromBeamspotJpsi.perp() / sqrt(jerr.rerr(displacementFromBeamspotJpsi)), PrescaleWeight);
0902               DiMuPVcos_.numerator->Fill(jpsi_cos, PrescaleWeight);
0903               DiMuPt_.numerator->Fill((m1.p4() + m.p4()).Pt(), PrescaleWeight);
0904               DiMuEta_.numerator->Fill((m1.p4() + m.p4()).Eta(), PrescaleWeight);
0905               DiMuPhi_.numerator->Fill((m1.p4() + m.p4()).Phi(), PrescaleWeight);
0906               DiMuDCA_.numerator->Fill(cApp.distance(), PrescaleWeight);
0907             }
0908 
0909             break;
0910 
0911           case 9:
0912 
0913             if (dimuonCL < minprob)
0914               continue;
0915             if (fabs(jpsi_cos) < mincos)
0916               continue;
0917             if ((displacementFromBeamspotJpsi.perp() / sqrt(jerr.rerr(displacementFromBeamspotJpsi))) < minDS)
0918               continue;
0919 
0920             if (trHandle.isValid())
0921               for (auto const& t : *trHandle) {
0922                 if (!trSelection_ref(t))
0923                   continue;
0924                 const reco::Track& itrk1 = t;
0925                 if (reco::deltaR(t, m1) <= min_dR)
0926                   continue;
0927                 if (reco::deltaR(t, m) <= min_dR)
0928                   continue;
0929                 if (!itrk1.quality(reco::TrackBase::highPurity))
0930                   continue;
0931 
0932                 // reconstruct B+ hadron
0933                 reco::Particle::LorentzVector pB, p1, p2, p3;
0934                 double trackMass2 = kaon_mass * kaon_mass;
0935                 double MuMass2 = mu_mass * mu_mass;
0936                 double e1 = sqrt(m.momentum().Mag2() + MuMass2);
0937                 double e2 = sqrt(m1.momentum().Mag2() + MuMass2);
0938                 double e3 = sqrt(itrk1.momentum().Mag2() + trackMass2);
0939                 p1 = reco::Particle::LorentzVector(m.px(), m.py(), m.pz(), e1);
0940                 p2 = reco::Particle::LorentzVector(m1.px(), m1.py(), m1.pz(), e2);
0941                 p3 = reco::Particle::LorentzVector(itrk1.px(), itrk1.py(), itrk1.pz(), e3);
0942                 pB = p1 + p2 + p3;
0943                 if (pB.mass() > maxmassJpsiTk || pB.mass() < minmassJpsiTk)
0944                   continue;
0945                 reco::TransientTrack trTT(itrk1, &magneticField);
0946                 std::vector<reco::TransientTrack> t_tks;
0947                 t_tks.push_back(mu1TT);
0948                 t_tks.push_back(mu2TT);
0949                 t_tks.push_back(trTT);
0950                 KalmanVertexFitter kvf;
0951                 TransientVertex tv = kvf.vertex(t_tks);
0952                 reco::Vertex vertex = tv;
0953                 if (!tv.isValid())
0954                   continue;
0955                 float JpsiTkCL = 0;
0956                 if ((vertex.chi2() >= 0.0) && (vertex.ndof() > 0))
0957                   JpsiTkCL = TMath::Prob(vertex.chi2(), vertex.ndof());
0958                 math::XYZVector pperp(m.px() + m1.px() + itrk1.px(), m.py() + m1.py() + itrk1.py(), 0.);
0959                 GlobalPoint secondaryVertex = tv.position();
0960                 GlobalError err = tv.positionError();
0961                 GlobalPoint displacementFromBeamspot(
0962                     -1 * ((vertexBeamSpot.x0() - secondaryVertex.x()) +
0963                           (secondaryVertex.z() - vertexBeamSpot.z0()) * vertexBeamSpot.dxdz()),
0964                     -1 * ((vertexBeamSpot.y0() - secondaryVertex.y()) +
0965                           (secondaryVertex.z() - vertexBeamSpot.z0()) * vertexBeamSpot.dydz()),
0966                     0);
0967                 reco::Vertex::Point vperp(displacementFromBeamspot.x(), displacementFromBeamspot.y(), 0.);
0968                 float jpsiKcos = vperp.Dot(pperp) / (vperp.R() * pperp.R());
0969                 if (JpsiTkCL < minprob)
0970                   continue;
0971                 if (fabs(jpsiKcos) < mincos)
0972                   continue;
0973                 if ((displacementFromBeamspot.perp() / sqrt(err.rerr(displacementFromBeamspot))) < minDS)
0974                   continue;
0975 
0976                 muPhi_.denominator->Fill(t.phi());
0977                 muEta_.denominator->Fill(t.eta());
0978                 muPt_.denominator->Fill(t.pt());
0979                 BMass_.denominator->Fill(pB.mass());
0980 
0981                 if (num_genTriggerEventFlag_->on() && num_genTriggerEventFlag_->accept(iEvent, iSetup)) {
0982                   if (!matchToTrigger(hltpath1, m1))
0983                     continue;
0984                   if (!matchToTrigger(hltpath1, m))
0985                     continue;
0986                   if (!matchToTrigger(hltpath1, t))
0987                     continue;
0988                   muPhi_.numerator->Fill(t.phi(), PrescaleWeight);
0989                   muEta_.numerator->Fill(t.eta(), PrescaleWeight);
0990                   muPt_.numerator->Fill(t.pt(), PrescaleWeight);
0991                   BMass_.numerator->Fill(pB.mass(), PrescaleWeight);
0992                 }
0993               }
0994 
0995             break;
0996 
0997           case 10:
0998 
0999             if (trHandle.isValid())
1000               for (auto const& t : *trHandle) {
1001                 if (!trSelection_ref(t))
1002                   continue;
1003                 const reco::Track& itrk1 = t;
1004                 if (reco::deltaR(t, m1) <= min_dR)
1005                   continue;
1006                 if (reco::deltaR(t, m) <= min_dR)
1007                   continue;
1008                 if (!itrk1.quality(reco::TrackBase::highPurity))
1009                   continue;
1010 
1011                 // reconstruct Mu+TkMu structure
1012                 reco::Particle::LorentzVector pB, p2, p3;
1013                 double trackMass2 = kaon_mass * kaon_mass;
1014                 double MuMass2 = mu_mass * mu_mass;
1015                 double e2 = sqrt(m1.momentum().Mag2() + MuMass2);
1016                 double e3 = sqrt(itrk1.momentum().Mag2() + trackMass2);
1017                 p2 = reco::Particle::LorentzVector(m1.px(), m1.py(), m1.pz(), e2);
1018                 p3 = reco::Particle::LorentzVector(itrk1.px(), itrk1.py(), itrk1.pz(), e3);
1019                 pB = p2 + p3;
1020                 if (pB.mass() > maxmassJpsiTk || pB.mass() < minmassJpsiTk)
1021                   continue;
1022                 reco::TransientTrack trTT(itrk1, &magneticField);
1023                 std::vector<reco::TransientTrack> t_tks;
1024                 t_tks.push_back(mu2TT);
1025                 t_tks.push_back(trTT);
1026                 KalmanVertexFitter kvf;
1027                 TransientVertex tv = kvf.vertex(t_tks);
1028                 reco::Vertex vertex = tv;
1029                 if (!tv.isValid())
1030                   continue;
1031                 float JpsiTkCL = 0;
1032                 if ((vertex.chi2() >= 0.0) && (vertex.ndof() > 0))
1033                   JpsiTkCL = TMath::Prob(vertex.chi2(), vertex.ndof());
1034                 math::XYZVector pperp(m1.px() + itrk1.px(), m1.py() + itrk1.py(), 0.);
1035                 GlobalPoint secondaryVertex = tv.position();
1036                 GlobalError err = tv.positionError();
1037                 GlobalPoint displacementFromBeamspot(
1038                     -1 * ((vertexBeamSpot.x0() - secondaryVertex.x()) +
1039                           (secondaryVertex.z() - vertexBeamSpot.z0()) * vertexBeamSpot.dxdz()),
1040                     -1 * ((vertexBeamSpot.y0() - secondaryVertex.y()) +
1041                           (secondaryVertex.z() - vertexBeamSpot.z0()) * vertexBeamSpot.dydz()),
1042                     0);
1043                 reco::Vertex::Point vperp(displacementFromBeamspot.x(), displacementFromBeamspot.y(), 0.);
1044                 if (JpsiTkCL < minprob)
1045                   continue;
1046 
1047                 muPhi_.denominator->Fill(m1.phi());
1048                 muEta_.denominator->Fill(m1.eta());
1049                 muPt_.denominator->Fill(m1.pt());
1050 
1051                 if (num_genTriggerEventFlag_->on() && num_genTriggerEventFlag_->accept(iEvent, iSetup)) {
1052                   if (!matchToTrigger(hltpath1, m1))
1053                     continue;
1054                   if (!matchToTrigger(hltpath1, m))
1055                     continue;
1056                   if (!matchToTrigger(hltpath1, t))
1057                     continue;
1058                   muPhi_.numerator->Fill(m1.phi(), PrescaleWeight);
1059                   muEta_.numerator->Fill(m1.eta(), PrescaleWeight);
1060                   muPt_.numerator->Fill(m1.pt(), PrescaleWeight);
1061                 }
1062               }
1063 
1064             break;
1065 
1066           case 11:
1067 
1068             if (dimuonCL < minprob)
1069               continue;
1070             if (fabs(jpsi_cos) < mincos)
1071               continue;
1072             if (displacementFromBeamspotJpsi.perp() / sqrt(jerr.rerr(displacementFromBeamspotJpsi)) < minDS)
1073               continue;
1074 
1075             if (trHandle.isValid())
1076               for (auto const& t : *trHandle) {
1077                 if (!trSelection_ref(t))
1078                   continue;
1079                 if ((reco::deltaR(t, m) <= min_dR))
1080                   continue;
1081                 if ((reco::deltaR(t, m1) <= min_dR))
1082                   continue;
1083 
1084                 for (auto const& t1 : *trHandle) {
1085                   if (&t - &(*trHandle)[0] >= &t1 - &(*trHandle)[0])
1086                     continue;  // not enough, need the following DeltaR checks
1087                   if (!trSelection_ref(t1))
1088                     continue;
1089                   if ((reco::deltaR(t1, m) <= min_dR))
1090                     continue;
1091                   if ((reco::deltaR(t1, m1) <= min_dR))
1092                     continue;
1093                   if ((reco::deltaR(t, t1) <= min_dR))
1094                     continue;
1095                   const reco::Track& itrk1 = t;
1096                   const reco::Track& itrk2 = t1;
1097                   if (!itrk1.quality(reco::TrackBase::highPurity))
1098                     continue;
1099                   if (!itrk2.quality(reco::TrackBase::highPurity))
1100                     continue;
1101 
1102                   // reconstruct Bs candidate
1103                   reco::Particle::LorentzVector pB, pTkTk, p1, p2, p3, p4;
1104                   double trackMass2 = kaon_mass * kaon_mass;
1105                   double MuMass2 = mu_mass * mu_mass;
1106                   double e1 = sqrt(m.momentum().Mag2() + MuMass2);
1107                   double e2 = sqrt(m1.momentum().Mag2() + MuMass2);
1108                   double e3 = sqrt(itrk1.momentum().Mag2() + trackMass2);
1109                   double e4 = sqrt(itrk2.momentum().Mag2() + trackMass2);
1110                   p1 = reco::Particle::LorentzVector(m.px(), m.py(), m.pz(), e1);
1111                   p2 = reco::Particle::LorentzVector(m1.px(), m1.py(), m1.pz(), e2);
1112                   p3 = reco::Particle::LorentzVector(itrk1.px(), itrk1.py(), itrk1.pz(), e3);
1113                   p4 = reco::Particle::LorentzVector(itrk2.px(), itrk2.py(), itrk2.pz(), e4);
1114                   pTkTk = p3 + p4;
1115                   if (pTkTk.mass() > maxmassTkTk || pTkTk.mass() < minmassTkTk)
1116                     continue;
1117                   pB = p1 + p2 + p3 + p4;
1118                   if (pB.mass() > maxmassJpsiTk || pB.mass() < minmassJpsiTk)
1119                     continue;
1120                   reco::TransientTrack mu1TT(m.track(), &magneticField);
1121                   reco::TransientTrack mu2TT(m1.track(), &magneticField);
1122                   reco::TransientTrack trTT(itrk1, &magneticField);
1123                   reco::TransientTrack tr1TT(itrk2, &magneticField);
1124                   std::vector<reco::TransientTrack> t_tks;
1125                   t_tks.push_back(mu1TT);
1126                   t_tks.push_back(mu2TT);
1127                   t_tks.push_back(trTT);
1128                   t_tks.push_back(tr1TT);
1129                   KalmanVertexFitter kvf;
1130                   TransientVertex tv = kvf.vertex(t_tks);  // this will compare the tracks
1131                   reco::Vertex vertex = tv;
1132                   if (!tv.isValid())
1133                     continue;
1134                   float JpsiTkCL = 0;
1135                   if ((vertex.chi2() >= 0.0) && (vertex.ndof() > 0))
1136                     JpsiTkCL = TMath::Prob(vertex.chi2(), vertex.ndof());
1137                   math::XYZVector pperp(
1138                       m.px() + m1.px() + itrk1.px() + itrk2.px(), m.py() + m1.py() + itrk1.py() + itrk2.py(), 0.);
1139                   GlobalPoint secondaryVertex = tv.position();
1140                   GlobalError err = tv.positionError();
1141                   GlobalPoint displacementFromBeamspot(
1142                       -1 * ((vertexBeamSpot.x0() - secondaryVertex.x()) +
1143                             (secondaryVertex.z() - vertexBeamSpot.z0()) * vertexBeamSpot.dxdz()),
1144                       -1 * ((vertexBeamSpot.y0() - secondaryVertex.y()) +
1145                             (secondaryVertex.z() - vertexBeamSpot.z0()) * vertexBeamSpot.dydz()),
1146                       0);
1147                   reco::Vertex::Point vperp(displacementFromBeamspot.x(), displacementFromBeamspot.y(), 0.);
1148                   float jpsiKcos = vperp.Dot(pperp) / (vperp.R() * pperp.R());
1149                   if (JpsiTkCL < minprob)
1150                     continue;
1151                   if (fabs(jpsiKcos) < mincos)
1152                     continue;
1153                   if ((displacementFromBeamspot.perp() / sqrt(err.rerr(displacementFromBeamspot))) < minDS)
1154                     continue;
1155 
1156                   mu1Phi_.denominator->Fill(t.phi());
1157                   mu1Eta_.denominator->Fill(t.eta());
1158                   mu1Pt_.denominator->Fill(t.pt());
1159                   mu2Phi_.denominator->Fill(t1.phi());
1160                   mu2Eta_.denominator->Fill(t1.eta());
1161                   mu2Pt_.denominator->Fill(t1.pt());
1162                   BMass_.denominator->Fill(pB.mass());
1163 
1164                   if (num_genTriggerEventFlag_->on() && num_genTriggerEventFlag_->accept(iEvent, iSetup)) {
1165                     if (!matchToTrigger(hltpath1, m))
1166                       continue;
1167                     if (!matchToTrigger(hltpath1, m1))
1168                       continue;
1169                     if (!matchToTrigger(hltpath1, t))
1170                       continue;
1171                     if (!matchToTrigger(hltpath1, t1))
1172                       continue;
1173                     mu1Phi_.numerator->Fill(t.phi(), PrescaleWeight);
1174                     mu1Eta_.numerator->Fill(t.eta(), PrescaleWeight);
1175                     mu1Pt_.numerator->Fill(t.pt(), PrescaleWeight);
1176                     mu2Phi_.numerator->Fill(t1.phi(), PrescaleWeight);
1177                     mu2Eta_.numerator->Fill(t1.eta(), PrescaleWeight);
1178                     mu2Pt_.numerator->Fill(t1.pt(), PrescaleWeight);
1179                     BMass_.numerator->Fill(pB.mass(), PrescaleWeight);
1180                   }
1181 
1182                 }  // for (auto const & t1 : *trHandle)
1183               }  // for (auto const & t : *trHandle)
1184 
1185             break;
1186         }
1187       }
1188     }
1189   }
1190 }
1191 
1192 void BPHMonitor::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
1193   edm::ParameterSetDescription desc;
1194   desc.add<std::string>("FolderName", "HLT/BPH/");
1195   desc.add<bool>("requireValidHLTPaths", true);
1196 
1197   desc.add<edm::InputTag>("tracks", edm::InputTag("generalTracks"));
1198   desc.add<edm::InputTag>("photons", edm::InputTag("photons"));
1199   desc.add<edm::InputTag>("offlinePVs", edm::InputTag("offlinePrimaryVertices"));
1200   desc.add<edm::InputTag>("beamSpot", edm::InputTag("offlineBeamSpot"));
1201   desc.add<edm::InputTag>("muons", edm::InputTag("muons"));
1202   desc.add<edm::InputTag>("hltTriggerSummaryAOD", edm::InputTag("hltTriggerSummaryAOD", "", "HLT"));
1203   desc.add<std::string>("muoSelection", "");
1204   desc.add<std::string>("muoSelection_ref",
1205                         "isPFMuon & isGlobalMuon  & innerTrack.hitPattern.trackerLayersWithMeasurement>5 & "
1206                         "innerTrack.hitPattern.numberOfValidPixelHits> 0");
1207   desc.add<std::string>(
1208       "muoSelection_tag",
1209       "isGlobalMuon && isPFMuon && isTrackerMuon && abs(eta) < 2.4 && innerTrack.hitPattern.numberOfValidPixelHits > 0 "
1210       "&& innerTrack.hitPattern.trackerLayersWithMeasurement > 5 && globalTrack.hitPattern.numberOfValidMuonHits > 0 "
1211       "&& globalTrack.normalizedChi2 < 10");  // tight selection for tag muon
1212   desc.add<std::string>("muoSelection_probe",
1213                         "isPFMuon & isGlobalMuon  & innerTrack.hitPattern.trackerLayersWithMeasurement>5 & "
1214                         "innerTrack.hitPattern.numberOfValidPixelHits> 0");
1215   desc.add<std::string>("trSelection_ref", "");
1216   desc.add<std::string>("DMSelection_ref", "Pt>4 & abs(eta)");
1217 
1218   desc.add<int>("nmuons", 1);
1219   desc.add<bool>("tnp", false);
1220   desc.add<int>("L3", 0);
1221   desc.add<int>("ptCut", 0);
1222   desc.add<int>("displaced", 0);
1223   desc.add<int>("trOrMu", 0);  // if =0, track param monitoring
1224   desc.add<int>("Jpsi", 0);
1225   desc.add<int>("Upsilon", 0);
1226   desc.add<int>("enum", 1);  // 1...9, 9 sets of variables to be filled, depends on the hlt path
1227   desc.add<int>("seagull", 1);
1228   desc.add<double>("maxmass", 3.596);
1229   desc.add<double>("minmass", 2.596);
1230   desc.add<double>("maxmassJpsi", 3.2);
1231   desc.add<double>("minmassJpsi", 3.);
1232   desc.add<double>("maxmassUpsilon", 10.0);
1233   desc.add<double>("minmassUpsilon", 8.8);
1234   desc.add<double>("maxmassTkTk", 10);
1235   desc.add<double>("minmassTkTk", 0);
1236   desc.add<double>("maxmassJpsiTk", 5.46);
1237   desc.add<double>("minmassJpsiTk", 5.1);
1238   desc.add<double>("kaon_mass", 0.493677);
1239   desc.add<double>("mu_mass", 0.1056583745);
1240   desc.add<double>("min_dR", 0.001);
1241   desc.add<double>("max_dR", 1.4);
1242   desc.add<double>("minprob", 0.005);
1243   desc.add<double>("mincos", 0.95);
1244   desc.add<double>("minDS", 3.);
1245   HLTPrescaleProvider::fillPSetDescription(
1246       desc, 1, edm::InputTag("gtStage2Digis"), edm::InputTag("gtStage2Digis"), false);
1247 
1248   edm::ParameterSetDescription genericTriggerEventPSet;
1249   GenericTriggerEventFlag::fillPSetDescription(genericTriggerEventPSet);
1250 
1251   desc.add<edm::ParameterSetDescription>("numGenericTriggerEventPSet", genericTriggerEventPSet);
1252   desc.add<edm::ParameterSetDescription>("denGenericTriggerEventPSet", genericTriggerEventPSet);
1253 
1254   edm::ParameterSetDescription histoPSet;
1255   edm::ParameterSetDescription phiPSet;
1256   edm::ParameterSetDescription etaPSet;
1257   edm::ParameterSetDescription ptPSet;
1258   edm::ParameterSetDescription dMu_ptPSet;
1259   edm::ParameterSetDescription d0PSet;
1260   edm::ParameterSetDescription z0PSet;
1261   edm::ParameterSetDescription dRPSet;
1262   edm::ParameterSetDescription massPSet;
1263   edm::ParameterSetDescription BmassPSet;
1264   edm::ParameterSetDescription dcaPSet;
1265   edm::ParameterSetDescription dsPSet;
1266   edm::ParameterSetDescription cosPSet;
1267   edm::ParameterSetDescription probPSet;
1268   edm::ParameterSetDescription TCoPSet;
1269   edm::ParameterSetDescription PUPSet;
1270   fillHistoPSetDescription(phiPSet);
1271   fillHistoPSetDescription(ptPSet);
1272   fillHistoPSetDescription(dMu_ptPSet);
1273   fillHistoPSetDescription(etaPSet);
1274   fillHistoPSetDescription(z0PSet);
1275   fillHistoPSetDescription(d0PSet);
1276   fillHistoPSetDescription(dRPSet);
1277   fillHistoPSetDescription(massPSet);
1278   fillHistoPSetDescription(BmassPSet);
1279   fillHistoPSetDescription(dcaPSet);
1280   fillHistoPSetDescription(dsPSet);
1281   fillHistoPSetDescription(cosPSet);
1282   fillHistoPSetDescription(probPSet);
1283   histoPSet.add<std::vector<double>>("ptBinning", {-0.5, 0, 2, 4, 8, 10, 12, 16, 20, 25, 30, 35, 40, 50});
1284   histoPSet.add<std::vector<double>>("dMuPtBinning", {6, 8, 12, 16, 20, 25, 30, 35, 40, 50, 70});
1285   histoPSet.add<std::vector<double>>("probBinning",
1286                                      {0.01, 0.02, 0.04, 0.06, 0.08, 0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1.0});
1287   histoPSet.add<edm::ParameterSetDescription>("d0PSet", d0PSet);
1288   histoPSet.add<edm::ParameterSetDescription>("etaPSet", etaPSet);
1289   histoPSet.add<edm::ParameterSetDescription>("phiPSet", phiPSet);
1290   histoPSet.add<edm::ParameterSetDescription>("z0PSet", z0PSet);
1291   histoPSet.add<edm::ParameterSetDescription>("dRPSet", dRPSet);
1292   histoPSet.add<edm::ParameterSetDescription>("massPSet", massPSet);
1293   histoPSet.add<edm::ParameterSetDescription>("BmassPSet", BmassPSet);
1294   histoPSet.add<edm::ParameterSetDescription>("dcaPSet", dcaPSet);
1295   histoPSet.add<edm::ParameterSetDescription>("dsPSet", dsPSet);
1296   histoPSet.add<edm::ParameterSetDescription>("cosPSet", cosPSet);
1297   desc.add<edm::ParameterSetDescription>("histoPSet", histoPSet);
1298 
1299   descriptions.add("bphMonitoring", desc);
1300 }
1301 
1302 std::string BPHMonitor::getTriggerName(std::string partialName) {
1303   const std::string trigger_name_tmp = partialName.substr(0, partialName.find("v*"));
1304   const unsigned int Ntriggers(hltConfig_.size());
1305   std::string trigger_name = "";
1306   for (unsigned int i = 0; i < Ntriggers; i++) {
1307     trigger_name = hltConfig_.triggerName(i);
1308     if (trigger_name.find(trigger_name_tmp) != std::string::npos)
1309       break;
1310   }
1311 
1312   return trigger_name;
1313 }
1314 
1315 template <typename T>
1316 bool BPHMonitor::matchToTrigger(const std::string& theTriggerName, T t) {
1317   bool matched = false;
1318   //validity check
1319   if (!hltConfig_.inited())
1320     return false;
1321 
1322   //Find the precise trigger name
1323   std::string trigger_name = getTriggerName(theTriggerName);
1324   const unsigned int trigger_index = hltConfig_.triggerIndex(trigger_name);
1325 
1326   //loop over all the modules for this trigger
1327   //by default use the last one
1328   unsigned int Nmodules = hltConfig_.size(trigger_index);
1329   const vector<string>& moduleLabels(hltConfig_.moduleLabels(trigger_index));
1330   unsigned int fIdx = 0;
1331   for (unsigned int i = 0; i < Nmodules; i++) {
1332     const unsigned int tmp_fIdx =
1333         handleTriggerEvent->filterIndex(edm::InputTag(moduleLabels[i], "", hltInputTag_1.process()));
1334     if (tmp_fIdx < handleTriggerEvent->sizeFilters())  //index of not used filters are set to sizeFilters()
1335     {
1336       fIdx = tmp_fIdx;
1337     }  //good index
1338   }
1339 
1340   //loop over all the objects in the filter of choice
1341   const trigger::Keys& KEYS(handleTriggerEvent->filterKeys(fIdx));
1342   const trigger::size_type nK(KEYS.size());
1343   const trigger::TriggerObjectCollection& TOC(handleTriggerEvent->getObjects());
1344   for (trigger::size_type i = 0; i != nK; ++i) {
1345     const trigger::TriggerObject& TO(TOC[KEYS[i]]);
1346     //perform matching: deltaR and pt check
1347     if ((reco::deltaR(t.eta(), t.phi(), TO.eta(), TO.phi()) <= 0.2) && (TMath::Abs(t.pt() - TO.pt()) < 0.12)) {
1348       matched = true;
1349     }
1350   }
1351 
1352   return matched;
1353 }
1354 
1355 double BPHMonitor::Prescale(const std::string hltpath1,
1356                             const std::string hltpath,
1357                             edm::Event const& iEvent,
1358                             edm::EventSetup const& iSetup,
1359                             HLTPrescaleProvider* hltPrescale_) {
1360   double Prescale_num = 1;
1361   double L1P = 1, HLTP = 1;
1362   bool flag = true;
1363   std::vector<bool> theSame_den;
1364   std::vector<bool> theSame_num;
1365   // object holding L1T (double) and HLT (uint) prescales
1366   auto const prescales_den = hltPrescale_->prescaleValuesInDetail<double>(iEvent, iSetup, hltpath);
1367   auto const prescales_num = hltPrescale_->prescaleValuesInDetail<double>(iEvent, iSetup, hltpath1);
1368   // retrieving HLT prescale
1369   auto PrescaleHLT_den = prescales_den.second;
1370   auto PrescaleHLT_num = prescales_num.second;
1371   if (PrescaleHLT_den > 0 && PrescaleHLT_num > 0)
1372     HLTP = PrescaleHLT_num / std::min(PrescaleHLT_num, PrescaleHLT_den);
1373 
1374   //retrieving L1 prescale
1375   //Checking if we have the same l1 seeds in den and num
1376   //taking into account that they can be written in different order in num and den
1377   //and some of them can be also switched off
1378 
1379   //check if for each den l1 there is the same l1 seed in num
1380   if (not prescales_den.first.empty()) {
1381     for (size_t iSeed = 0; iSeed < prescales_den.first.size(); ++iSeed) {
1382       auto l1_den = prescales_den.first.at(iSeed).first;
1383       auto l1_denp = prescales_den.first.at(iSeed).second;
1384       if (l1_denp < 1)
1385         continue;
1386       flag = false;
1387       for (size_t iSeed1 = 0; iSeed1 < prescales_num.first.size(); ++iSeed1) {
1388         auto l1_num = prescales_num.first.at(iSeed1).first;
1389         auto l1_nump = prescales_num.first.at(iSeed1).second;
1390         if (l1_num == l1_den && l1_nump >= 1)  //the same seed
1391         {
1392           flag = true;
1393           break;
1394         }
1395       }
1396       theSame_den.push_back(flag);
1397     }
1398   }
1399   //check if for each num l1 there is the same l1 seed in den
1400   if (not prescales_num.first.empty()) {
1401     for (size_t iSeed = 0; iSeed < prescales_num.first.size(); ++iSeed) {
1402       auto l1_num = prescales_num.first.at(iSeed).first;
1403       auto l1_nump = prescales_num.first.at(iSeed).second;
1404       if (l1_nump < 1)
1405         continue;
1406       flag = false;
1407       for (size_t iSeed1 = 0; iSeed1 < prescales_den.first.size(); ++iSeed1) {
1408         auto l1_den = prescales_den.first.at(iSeed1).first;
1409         auto l1_denp = prescales_den.first.at(iSeed1).second;
1410         if (l1_den == l1_num && l1_denp >= 1)  //the same seed
1411         {
1412           flag = true;
1413           break;
1414         }
1415       }
1416       theSame_num.push_back(flag);
1417     }
1418   }
1419   flag = true;
1420 
1421   if (theSame_num.size() == theSame_den.size()) {
1422     for (size_t i = 0; i < theSame_num.size(); ++i) {
1423       if ((!theSame_num.at(i)) || (!theSame_den.at(i))) {
1424         flag = false;
1425         break;
1426       }
1427     }
1428   }
1429 
1430   if (flag && (theSame_num.size() == theSame_den.size())) {
1431     L1P = 1;  //den and num have the same set of l1 seeds
1432   } else {
1433     if (not prescales_num.first.empty()) {
1434       Prescale_num = 1;
1435       for (size_t iSeed = 0; iSeed < prescales_num.first.size(); ++iSeed) {
1436         auto l1 = prescales_num.first.at(iSeed).second;
1437         if (l1 < 1)
1438           continue;
1439         if (l1 == 1) {
1440           Prescale_num = 1;
1441           break;
1442         } else
1443           Prescale_num *= 1 - (1.0 / (l1));
1444       }
1445       if (Prescale_num != 1)
1446         Prescale_num = 1.0 / (1 - Prescale_num);
1447     }
1448     L1P = Prescale_num;
1449   }
1450 
1451   return L1P * HLTP;
1452 }
1453 
1454 DEFINE_FWK_MODULE(BPHMonitor);