Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:59

0001 // -----------------------------
0002 //
0003 // Offline DQM for razor triggers. The razor inclusive analysis measures trigger efficiency
0004 // in SingleElectron events (orthogonal to analysis), as a 2D function of the razor variables
0005 // M_R and R^2. Also monitor dPhi_R, used offline for  QCD and/or detector-related MET tail
0006 // rejection.
0007 // Based on DQMOffline/Trigger/plugins/METMonitor.cc
0008 //
0009 // -----------------------------
0010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0011 #include "FWCore/Framework/interface/Frameworkfwd.h"
0012 #include "FWCore/Framework/interface/Event.h"
0013 #include "FWCore/Framework/interface/MakerMacros.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0015 #include "FWCore/ServiceRegistry/interface/Service.h"
0016 #include "DQMServices/Core/interface/DQMStore.h"
0017 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0018 #include "DQMOffline/Trigger/plugins/TriggerDQMBase.h"
0019 #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
0020 #include "CommonTools/TriggerUtils/interface/GenericTriggerEventFlag.h"
0021 #include "DataFormats/METReco/interface/PFMET.h"
0022 #include "DataFormats/METReco/interface/PFMETCollection.h"
0023 #include "DataFormats/JetReco/interface/PFJet.h"
0024 #include "DataFormats/JetReco/interface/PFJetCollection.h"
0025 #include "HLTrigger/JetMET/interface/HLTRHemisphere.h"
0026 #include "DataFormats/Math/interface/LorentzVector.h"
0027 #include "DataFormats/Math/interface/deltaPhi.h"
0028 
0029 #include <TVector3.h>
0030 
0031 class RazorMonitor : public DQMEDAnalyzer, public TriggerDQMBase {
0032 public:
0033   typedef dqm::reco::MonitorElement MonitorElement;
0034   typedef dqm::reco::DQMStore DQMStore;
0035 
0036   RazorMonitor(const edm::ParameterSet&);
0037   ~RazorMonitor() throw() override;
0038   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0039 
0040   static double CalcMR(const math::XYZTLorentzVector& ja, const math::XYZTLorentzVector& jb);
0041   static double CalcR(double MR,
0042                       const math::XYZTLorentzVector& ja,
0043                       const math::XYZTLorentzVector& jb,
0044                       const edm::Handle<std::vector<reco::PFMET> >& met);
0045 
0046 protected:
0047   void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
0048   void analyze(edm::Event const& iEvent, edm::EventSetup const& iSetup) override;
0049 
0050 private:
0051   const std::string folderName_;
0052 
0053   const bool requireValidHLTPaths_;
0054   bool hltPathsAreValid_;
0055 
0056   edm::EDGetTokenT<reco::PFMETCollection> metToken_;
0057   edm::EDGetTokenT<reco::PFJetCollection> jetToken_;
0058   edm::EDGetTokenT<std::vector<math::XYZTLorentzVector> > theHemispheres_;
0059 
0060   std::vector<double> rsq_binning_;
0061   std::vector<double> mr_binning_;
0062   std::vector<double> dphiR_binning_;
0063 
0064   ObjME MR_ME_;
0065   ObjME Rsq_ME_;
0066   ObjME dPhiR_ME_;
0067   ObjME MRVsRsq_ME_;
0068 
0069   std::unique_ptr<GenericTriggerEventFlag> num_genTriggerEventFlag_;
0070   std::unique_ptr<GenericTriggerEventFlag> den_genTriggerEventFlag_;
0071 
0072   StringCutObjectSelector<reco::MET, true> metSelection_;
0073   StringCutObjectSelector<reco::PFJet, true> jetSelection_;
0074   unsigned int njets_;
0075   float rsqCut_;
0076   float mrCut_;
0077 };
0078 
0079 // -----------------------------
0080 //
0081 // Offline DQM for razor triggers. The razor inclusive analysis measures trigger efficiency
0082 // in SingleElectron events (orthogonal to analysis), as a 2D function of the razor variables
0083 // M_R and R^2. Also monitor dPhi_R, used offline for  QCD and/or detector-related MET tail
0084 // rejection.
0085 // Based on DQMOffline/Trigger/plugins/METMonitor.*
0086 //
0087 // -----------------------------
0088 
0089 RazorMonitor::RazorMonitor(const edm::ParameterSet& iConfig)
0090     : folderName_(iConfig.getParameter<std::string>("FolderName")),
0091       requireValidHLTPaths_(iConfig.getParameter<bool>("requireValidHLTPaths")),
0092       hltPathsAreValid_(false),
0093       metToken_(consumes<reco::PFMETCollection>(iConfig.getParameter<edm::InputTag>("met"))),
0094       jetToken_(mayConsume<reco::PFJetCollection>(iConfig.getParameter<edm::InputTag>("jets"))),
0095       theHemispheres_(
0096           consumes<std::vector<math::XYZTLorentzVector> >(iConfig.getParameter<edm::InputTag>("hemispheres"))),
0097       rsq_binning_(iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("rsqBins")),
0098       mr_binning_(iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("mrBins")),
0099       dphiR_binning_(
0100           iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("dphiRBins")),
0101       num_genTriggerEventFlag_(new GenericTriggerEventFlag(
0102           iConfig.getParameter<edm::ParameterSet>("numGenericTriggerEventPSet"), consumesCollector(), *this)),
0103       den_genTriggerEventFlag_(new GenericTriggerEventFlag(
0104           iConfig.getParameter<edm::ParameterSet>("denGenericTriggerEventPSet"), consumesCollector(), *this)),
0105       metSelection_(iConfig.getParameter<std::string>("metSelection")),
0106       jetSelection_(iConfig.getParameter<std::string>("jetSelection")),
0107       njets_(iConfig.getParameter<unsigned int>("njets")),
0108       rsqCut_(iConfig.getParameter<double>("rsqCut")),
0109       mrCut_(iConfig.getParameter<double>("mrCut")) {}
0110 
0111 RazorMonitor::~RazorMonitor() throw() {
0112   if (num_genTriggerEventFlag_) {
0113     num_genTriggerEventFlag_.reset();
0114   }
0115   if (den_genTriggerEventFlag_) {
0116     den_genTriggerEventFlag_.reset();
0117   }
0118 }
0119 
0120 void RazorMonitor::bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& iRun, edm::EventSetup const& iSetup) {
0121   // Initialize the GenericTriggerEventFlag
0122   if (num_genTriggerEventFlag_ && num_genTriggerEventFlag_->on())
0123     num_genTriggerEventFlag_->initRun(iRun, iSetup);
0124   if (den_genTriggerEventFlag_ && den_genTriggerEventFlag_->on())
0125     den_genTriggerEventFlag_->initRun(iRun, iSetup);
0126 
0127   // check if every HLT path specified in numerator and denominator has a valid match in the HLT Menu
0128   hltPathsAreValid_ = (num_genTriggerEventFlag_ && den_genTriggerEventFlag_ && num_genTriggerEventFlag_->on() &&
0129                        den_genTriggerEventFlag_->on() && num_genTriggerEventFlag_->allHLTPathsAreValid() &&
0130                        den_genTriggerEventFlag_->allHLTPathsAreValid());
0131 
0132   // if valid HLT paths are required,
0133   // create DQM outputs only if all paths are valid
0134   if (requireValidHLTPaths_ and (not hltPathsAreValid_)) {
0135     return;
0136   }
0137 
0138   std::string histname, histtitle;
0139 
0140   std::string currentFolder = folderName_;
0141   ibooker.setCurrentFolder(currentFolder);
0142 
0143   // 1D hist, MR
0144   histname = "MR";
0145   histtitle = "PF MR";
0146   bookME(ibooker, MR_ME_, histname, histtitle, mr_binning_);
0147   setMETitle(MR_ME_, "PF M_{R} [GeV]", "events / [GeV]");
0148 
0149   // 1D hist, Rsq
0150   histname = "Rsq";
0151   histtitle = "PF Rsq";
0152   bookME(ibooker, Rsq_ME_, histname, histtitle, rsq_binning_);
0153   setMETitle(Rsq_ME_, "PF R^{2}", "events");
0154 
0155   // 1D hist, dPhiR
0156   histname = "dPhiR";
0157   histtitle = "dPhiR";
0158   bookME(ibooker, dPhiR_ME_, histname, histtitle, dphiR_binning_);
0159   setMETitle(dPhiR_ME_, "dPhi_{R}", "events");
0160 
0161   // 2D hist, MR & Rsq
0162   histname = "MRVsRsq";
0163   histtitle = "PF MR vs PF Rsq";
0164   bookME(ibooker, MRVsRsq_ME_, histname, histtitle, mr_binning_, rsq_binning_);
0165   setMETitle(MRVsRsq_ME_, "M_{R} [GeV]", "R^{2}");
0166 }
0167 
0168 void RazorMonitor::analyze(edm::Event const& iEvent, edm::EventSetup const& iSetup) {
0169   // if valid HLT paths are required,
0170   // analyze event only if all paths are valid
0171   if (requireValidHLTPaths_ and (not hltPathsAreValid_)) {
0172     return;
0173   }
0174 
0175   // Filter out events if Trigger Filtering is requested
0176   if (den_genTriggerEventFlag_->on() && !den_genTriggerEventFlag_->accept(iEvent, iSetup))
0177     return;
0178 
0179   //met collection
0180   edm::Handle<reco::PFMETCollection> metHandle;
0181   iEvent.getByToken(metToken_, metHandle);
0182   reco::PFMET pfmet = metHandle->front();
0183   if (!metSelection_(pfmet))
0184     return;
0185 
0186   //jet collection, track # of jets for two working points
0187   edm::Handle<reco::PFJetCollection> jetHandle;
0188   iEvent.getByToken(jetToken_, jetHandle);
0189   std::vector<reco::PFJet> jets;
0190   if (jetHandle->size() < njets_)
0191     return;
0192   for (auto const& j : *jetHandle) {
0193     if (jetSelection_(j))
0194       jets.push_back(j);
0195   }
0196   if (jets.size() < njets_)
0197     return;
0198 
0199   //razor hemisphere clustering from previous step
0200   edm::Handle<vector<math::XYZTLorentzVector> > hemispheres;
0201   iEvent.getByToken(theHemispheres_, hemispheres);
0202 
0203   if (not hemispheres.isValid()) {
0204     return;
0205   }
0206 
0207   if (hemispheres
0208           ->empty()) {  // the Hemisphere Maker will produce an empty collection of hemispheres if # of jets is too high
0209     edm::LogError("DQM_HLT_Razor") << "Cannot calculate M_R and R^2 because there are too many jets! (trigger passed "
0210                                       "automatically without forming the hemispheres)"
0211                                    << endl;
0212     return;
0213   }
0214 
0215   // should always have 2 hemispheres -- no muons included (c. 2017), if not return invalid hemisphere collection
0216   // retaining check for hemisphere size 5 or 10 which correspond to the one or two muon case for possible future use
0217   if (!hemispheres->empty() && hemispheres->size() != 2 && hemispheres->size() != 5 && hemispheres->size() != 10) {
0218     edm::LogError("DQM_HLT_Razor") << "Invalid hemisphere collection!  hemispheres->size() = " << hemispheres->size()
0219                                    << endl;
0220     return;
0221   }
0222 
0223   //calculate razor variables, with hemispheres pT-ordered
0224   double MR = 0, R = 0;
0225   if (hemispheres->at(1).Pt() > hemispheres->at(0).Pt()) {
0226     MR = CalcMR(hemispheres->at(1), hemispheres->at(0));
0227     R = CalcR(MR, hemispheres->at(1), hemispheres->at(0), metHandle);
0228   } else {
0229     MR = CalcMR(hemispheres->at(0), hemispheres->at(1));
0230     R = CalcR(MR, hemispheres->at(0), hemispheres->at(1), metHandle);
0231   }
0232 
0233   double Rsq = R * R;
0234   double dPhiR = abs(deltaPhi(hemispheres->at(0).Phi(), hemispheres->at(1).Phi()));
0235 
0236   //apply offline selection cuts
0237   if (Rsq < rsqCut_ && MR < mrCut_)
0238     return;
0239 
0240   // applying selection for denominator
0241   if (den_genTriggerEventFlag_->on() && !den_genTriggerEventFlag_->accept(iEvent, iSetup))
0242     return;
0243 
0244   // filling histograms (denominator)
0245   if (Rsq >= rsqCut_) {
0246     MR_ME_.denominator->Fill(MR);
0247   }
0248 
0249   if (MR >= mrCut_) {
0250     Rsq_ME_.denominator->Fill(Rsq);
0251   }
0252 
0253   dPhiR_ME_.denominator->Fill(dPhiR);
0254 
0255   MRVsRsq_ME_.denominator->Fill(MR, Rsq);
0256 
0257   // applying selection for numerator
0258   if (num_genTriggerEventFlag_->on() && !num_genTriggerEventFlag_->accept(iEvent, iSetup))
0259     return;
0260 
0261   // filling histograms (numerator)
0262   if (Rsq >= rsqCut_) {
0263     MR_ME_.numerator->Fill(MR);
0264   }
0265 
0266   if (MR >= mrCut_) {
0267     Rsq_ME_.numerator->Fill(Rsq);
0268   }
0269 
0270   dPhiR_ME_.numerator->Fill(dPhiR);
0271 
0272   MRVsRsq_ME_.numerator->Fill(MR, Rsq);
0273 }
0274 
0275 void RazorMonitor::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0276   edm::ParameterSetDescription desc;
0277   desc.add<std::string>("FolderName", "HLT/SUSY/Razor");
0278   desc.add<bool>("requireValidHLTPaths", true);
0279 
0280   desc.add<edm::InputTag>("met", edm::InputTag("pfMet"));
0281   desc.add<edm::InputTag>("jets", edm::InputTag("ak4PFJetsCHS"));
0282   desc.add<edm::InputTag>("hemispheres", edm::InputTag("hemispheresDQM"))
0283       ->setComment("hemisphere jets used to compute razor variables");
0284   desc.add<std::string>("metSelection", "pt > 0");
0285 
0286   // from 2016 offline selection
0287   desc.add<std::string>("jetSelection", "pt > 80");
0288   desc.add<unsigned int>("njets", 2);
0289   desc.add<double>("mrCut", 300);
0290   desc.add<double>("rsqCut", 0.15);
0291 
0292   edm::ParameterSetDescription genericTriggerEventPSet;
0293   GenericTriggerEventFlag::fillPSetDescription(genericTriggerEventPSet);
0294   desc.add<edm::ParameterSetDescription>("numGenericTriggerEventPSet", genericTriggerEventPSet);
0295   desc.add<edm::ParameterSetDescription>("denGenericTriggerEventPSet", genericTriggerEventPSet);
0296 
0297   //binning from 2016 offline selection
0298   edm::ParameterSetDescription histoPSet;
0299   std::vector<double> mrbins = {0., 100., 200., 300., 400., 500., 575., 650., 750., 900., 1200., 1600., 2500., 4000.};
0300   histoPSet.add<std::vector<double> >("mrBins", mrbins);
0301 
0302   std::vector<double> rsqbins = {0., 0.05, 0.1, 0.15, 0.2, 0.25, 0.30, 0.41, 0.52, 0.64, 0.8, 1.5};
0303   histoPSet.add<std::vector<double> >("rsqBins", rsqbins);
0304 
0305   std::vector<double> dphirbins = {0., 0.5, 1.0, 1.5, 2.0, 2.5, 2.8, 3.0, 3.2};
0306   histoPSet.add<std::vector<double> >("dphiRBins", dphirbins);
0307 
0308   desc.add<edm::ParameterSetDescription>("histoPSet", histoPSet);
0309 
0310   descriptions.add("razorMonitoring", desc);
0311 }
0312 
0313 //CalcMR and CalcR borrowed from HLTRFilter.cc
0314 double RazorMonitor::CalcMR(const math::XYZTLorentzVector& ja, const math::XYZTLorentzVector& jb) {
0315   if (ja.Pt() <= 0.1)
0316     return -1;
0317 
0318   double A = ja.P();
0319   double B = jb.P();
0320   double az = ja.Pz();
0321   double bz = jb.Pz();
0322   TVector3 jaT, jbT;
0323   jaT.SetXYZ(ja.Px(), ja.Py(), 0.0);
0324   jbT.SetXYZ(jb.Px(), jb.Py(), 0.0);
0325   double ATBT = (jaT + jbT).Mag2();
0326 
0327   double MR = sqrt((A + B) * (A + B) - (az + bz) * (az + bz) -
0328                    (jbT.Dot(jbT) - jaT.Dot(jaT)) * (jbT.Dot(jbT) - jaT.Dot(jaT)) / (jaT + jbT).Mag2());
0329   double mybeta = (jbT.Dot(jbT) - jaT.Dot(jaT)) / sqrt(ATBT * ((A + B) * (A + B) - (az + bz) * (az + bz)));
0330 
0331   double mygamma = 1. / sqrt(1. - mybeta * mybeta);
0332 
0333   //use gamma times MRstar
0334   return MR * mygamma;
0335 }
0336 
0337 double RazorMonitor::CalcR(double MR,
0338                            const math::XYZTLorentzVector& ja,
0339                            const math::XYZTLorentzVector& jb,
0340                            const edm::Handle<std::vector<reco::PFMET> >& inputMet) {
0341   //now we can calculate MTR
0342   const math::XYZVector met = (inputMet->front()).momentum();
0343 
0344   double MTR = sqrt(0.5 * (met.R() * (ja.Pt() + jb.Pt()) - met.Dot(ja.Vect() + jb.Vect())));
0345 
0346   //filter events
0347   return float(MTR) / float(MR);  //R
0348 }
0349 
0350 // Define this as a plug-in
0351 DEFINE_FWK_MODULE(RazorMonitor);