Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "DQM/DataScouting/plugins/DiJetVarAnalyzer.h"
0002 
0003 #include "DataFormats/EgammaCandidates/interface/Electron.h"
0004 #include "DataFormats/EgammaCandidates/interface/ElectronFwd.h"
0005 #include "DataFormats/JetReco/interface/CaloJet.h"
0006 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
0007 #include "DataFormats/METReco/interface/CaloMET.h"
0008 #include "DataFormats/METReco/interface/CaloMETCollection.h"
0009 #include "DataFormats/MuonReco/interface/Muon.h"
0010 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
0011 
0012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0013 
0014 #include "HLTrigger/HLTcore/interface/TriggerExpressionData.h"
0015 #include "HLTrigger/HLTcore/interface/TriggerExpressionEvaluator.h"
0016 #include "HLTrigger/HLTcore/interface/TriggerExpressionParser.h"
0017 
0018 #include <cmath>
0019 
0020 // A simple constructor which takes as inoput only the name of the PF jet
0021 // collection
0022 DiJetVarAnalyzer::DiJetVarAnalyzer(const edm::ParameterSet &conf)
0023     : ScoutingAnalyzerBase(conf),
0024       jetCollectionTag_(conf.getUntrackedParameter<edm::InputTag>("jetCollectionTag")),
0025       widejetsCollectionTag_(conf.getUntrackedParameter<edm::InputTag>("widejetsCollectionTag")),
0026       metCollectionTag_(conf.getUntrackedParameter<edm::InputTag>("metCollectionTag")),
0027       metCleanCollectionTag_(conf.getUntrackedParameter<edm::InputTag>("metCleanCollectionTag")),
0028       numwidejets_(conf.getParameter<unsigned int>("numwidejets")),
0029       etawidejets_(conf.getParameter<double>("etawidejets")),
0030       ptwidejets_(conf.getParameter<double>("ptwidejets")),
0031       detawidejets_(conf.getParameter<double>("detawidejets")),
0032       dphiwidejets_(conf.getParameter<double>("dphiwidejets")),
0033       maxEMfraction_(conf.getParameter<double>("maxEMfraction")),
0034       maxHADfraction_(conf.getParameter<double>("maxHADfraction")),
0035       HLTpathMain_(triggerExpression::parse(conf.getParameter<std::string>("HLTpathMain"))),
0036       HLTpathMonitor_(triggerExpression::parse(conf.getParameter<std::string>("HLTpathMonitor"))),
0037       triggerConfiguration_(conf.getParameterSet("triggerConfiguration"), consumesCollector()) {
0038   // set Token(-s)
0039   jetCollectionTagToken_ =
0040       consumes<reco::CaloJetCollection>(conf.getUntrackedParameter<edm::InputTag>("jetCollectionTag"));
0041   widejetsCollectionTagToken_ = consumes<std::vector<math::PtEtaPhiMLorentzVector>>(
0042       conf.getUntrackedParameter<edm::InputTag>("widejetsCollectionTag"));
0043   metCollectionTagToken_ =
0044       consumes<reco::CaloMETCollection>(conf.getUntrackedParameter<edm::InputTag>("metCollectionTag"));
0045   metCleanCollectionTagToken_ =
0046       consumes<reco::CaloMETCollection>(conf.getUntrackedParameter<edm::InputTag>("metCleanCollectionTag"));
0047 }
0048 
0049 DiJetVarAnalyzer::~DiJetVarAnalyzer() {}
0050 
0051 // Function to book the Monitoring Elements.
0052 void DiJetVarAnalyzer::bookHistograms(DQMStore::IBooker &iBooker, edm::Run const &, edm::EventSetup const &) {
0053   ScoutingAnalyzerBase::prepareBooking(iBooker);
0054   // ==> TO BE UPDATED FOR sqrt(s)=8 TeV
0055   const int N_mass_bins = 83;
0056   float massBins[N_mass_bins + 1] = {
0057       1,    3,    6,    10,   16,   23,   31,   40,   50,   61,   74,   88,   103,  119,  137,  156,  176,
0058       197,  220,  244,  270,  296,  325,  354,  386,  419,  453,  489,  526,  565,  606,  649,  693,  740,
0059       788,  838,  890,  944,  1000, 1058, 1118, 1181, 1246, 1313, 1383, 1455, 1530, 1607, 1687, 1770, 1856,
0060       1945, 2037, 2132, 2231, 2332, 2438, 2546, 2659, 2775, 2895, 3019, 3147, 3279, 3416, 3558, 3704, 3854,
0061       4010, 4171, 4337, 4509, 4686, 4869, 5058, 5253, 5455, 5663, 5877, 6099, 6328, 6564, 6808, 7000};
0062 
0063   // 1D histograms
0064   m_cutFlow = bookH1withSumw2(iBooker, "h1_cutFlow", "Cut Flow", 7, 0., 7., "Cut", "Number of events");
0065   m_cutFlow->getTH1()->GetXaxis()->SetBinLabel(1, "No cut");
0066   m_cutFlow->getTH1()->GetXaxis()->SetBinLabel(2, "N(WideJets)>=2");
0067   m_cutFlow->getTH1()->GetXaxis()->SetBinLabel(3, "|#eta|<2.5 , pT>30");
0068   m_cutFlow->getTH1()->GetXaxis()->SetBinLabel(4, "|#Delta#eta|<1.3");
0069   m_cutFlow->getTH1()->GetXaxis()->SetBinLabel(5, "JetID");
0070   m_cutFlow->getTH1()->GetXaxis()->SetBinLabel(6, "|#Delta#phi|>#pi/3");
0071   m_cutFlow->getTH1()->GetXaxis()->SetBinLabel(7, "|met-metClean|>0.1");
0072 
0073   m_MjjWide_finalSel = bookH1withSumw2(iBooker,
0074                                        "h1_MjjWide_finalSel",
0075                                        "M_{jj} WideJets (final selection)",
0076                                        8000,
0077                                        0.,
0078                                        8000.,
0079                                        "M_{jj} WideJets [GeV]",
0080                                        "Number of events");
0081 
0082   m_MjjWide_finalSel_varbin = bookH1withSumw2BinArray(iBooker,
0083                                                       "h1_MjjWide_finalSel_varbin",
0084                                                       "M_{jj} WideJets (final selection)",
0085                                                       N_mass_bins,
0086                                                       massBins,
0087                                                       "M_{jj} WideJets [GeV]",
0088                                                       "Number of events");
0089 
0090   m_MjjWide_finalSel_WithoutNoiseFilter = bookH1withSumw2(iBooker,
0091                                                           "h1_MjjWide_finalSel_WithoutNoiseFilter",
0092                                                           "M_{jj} WideJets (final selection, without noise filters)",
0093                                                           8000,
0094                                                           0.,
0095                                                           8000.,
0096                                                           "M_{jj} WideJets [GeV]",
0097                                                           "Number of events");
0098 
0099   m_MjjWide_finalSel_WithoutNoiseFilter_varbin =
0100       bookH1withSumw2BinArray(iBooker,
0101                               "h1_MjjWide_finalSel_WithoutNoiseFilter_varbin",
0102                               "M_{jj} WideJets (final selection, without noise filters)",
0103                               N_mass_bins,
0104                               massBins,
0105                               "M_{jj} WideJets [GeV]",
0106                               "Number of events");
0107 
0108   m_MjjWide_deta_0p0_0p5 = bookH1withSumw2(iBooker,
0109                                            "h1_MjjWide_deta_0p0_0p5",
0110                                            "M_{jj} WideJets (0.0<=#Delta#eta<0.5)",
0111                                            8000,
0112                                            0.,
0113                                            8000.,
0114                                            "M_{jj} WideJets [GeV]",
0115                                            "Number of events");
0116 
0117   m_MjjWide_deta_0p5_1p0 = bookH1withSumw2(iBooker,
0118                                            "h1_MjjWide_deta_0p5_1p0",
0119                                            "M_{jj} WideJets (0.5<=#Delta#eta<1.0)",
0120                                            8000,
0121                                            0.,
0122                                            8000.,
0123                                            "M_{jj} WideJets [GeV]",
0124                                            "Number of events");
0125 
0126   m_MjjWide_deta_1p0_1p5 = bookH1withSumw2(iBooker,
0127                                            "h1_MjjWide_deta_1p0_1p5",
0128                                            "M_{jj} WideJets (1.0<=#Delta#eta<1.5)",
0129                                            8000,
0130                                            0.,
0131                                            8000.,
0132                                            "M_{jj} WideJets [GeV]",
0133                                            "Number of events");
0134 
0135   m_MjjWide_deta_1p5_2p0 = bookH1withSumw2(iBooker,
0136                                            "h1_MjjWide_deta_1p5_2p0",
0137                                            "M_{jj} WideJets (1.5<=#Delta#eta<2.0)",
0138                                            8000,
0139                                            0.,
0140                                            8000.,
0141                                            "M_{jj} WideJets [GeV]",
0142                                            "Number of events");
0143 
0144   m_MjjWide_deta_2p0_2p5 = bookH1withSumw2(iBooker,
0145                                            "h1_MjjWide_deta_2p0_2p5",
0146                                            "M_{jj} WideJets (2.0<=#Delta#eta<2.5)",
0147                                            8000,
0148                                            0.,
0149                                            8000.,
0150                                            "M_{jj} WideJets [GeV]",
0151                                            "Number of events");
0152 
0153   m_MjjWide_deta_2p5_3p0 = bookH1withSumw2(iBooker,
0154                                            "h1_MjjWide_deta_2p5_3p0",
0155                                            "M_{jj} WideJets (2.5<=#Delta#eta<3.0)",
0156                                            8000,
0157                                            0.,
0158                                            8000.,
0159                                            "M_{jj} WideJets [GeV]",
0160                                            "Number of events");
0161 
0162   m_MjjWide_deta_3p0_inf = bookH1withSumw2(iBooker,
0163                                            "h1_MjjWide_deta_3p0_inf",
0164                                            "M_{jj} WideJets (#Delta#eta>=3.0)",
0165                                            8000,
0166                                            0.,
0167                                            8000.,
0168                                            "M_{jj} WideJets [GeV]",
0169                                            "Number of events");
0170 
0171   m_MjjWide_den_NOdeta = bookH1withSumw2(iBooker,
0172                                          "h1_MjjWide_den_NOdeta",
0173                                          "HLT Efficiency Studies (no deta cut)",
0174                                          400,
0175                                          0.,
0176                                          2000.,
0177                                          "M_{jj} WideJets [GeV]",
0178                                          "Number of events");
0179 
0180   m_MjjWide_num_NOdeta = bookH1withSumw2(iBooker,
0181                                          "h1_MjjWide_num_NOdeta",
0182                                          "HLT Efficiency Studies (no deta cut)",
0183                                          400,
0184                                          0.,
0185                                          2000.,
0186                                          "M_{jj} WideJets [GeV]",
0187                                          "Number of events");
0188 
0189   m_MjjWide_den_detaL4 = bookH1withSumw2(iBooker,
0190                                          "h1_MjjWide_den_detaL4",
0191                                          "HLT Efficiency Studies (deta cut < 4.0)",
0192                                          400,
0193                                          0.,
0194                                          2000.,
0195                                          "M_{jj} WideJets [GeV]",
0196                                          "Number of events");
0197 
0198   m_MjjWide_num_detaL4 = bookH1withSumw2(iBooker,
0199                                          "h1_MjjWide_num_detaL4",
0200                                          "HLT Efficiency Studies (deta cut < 4.0)",
0201                                          400,
0202                                          0.,
0203                                          2000.,
0204                                          "M_{jj} WideJets [GeV]",
0205                                          "Number of events");
0206 
0207   m_MjjWide_den_detaL3 = bookH1withSumw2(iBooker,
0208                                          "h1_MjjWide_den_detaL3",
0209                                          "HLT Efficiency Studies (deta cut < 3.0)",
0210                                          400,
0211                                          0.,
0212                                          2000.,
0213                                          "M_{jj} WideJets [GeV]",
0214                                          "Number of events");
0215 
0216   m_MjjWide_num_detaL3 = bookH1withSumw2(iBooker,
0217                                          "h1_MjjWide_num_detaL3",
0218                                          "HLT Efficiency Studies (deta cut < 3.0)",
0219                                          400,
0220                                          0.,
0221                                          2000.,
0222                                          "M_{jj} WideJets [GeV]",
0223                                          "Number of events");
0224 
0225   m_MjjWide_den_detaL2 = bookH1withSumw2(iBooker,
0226                                          "h1_MjjWide_den_detaL2",
0227                                          "HLT Efficiency Studies (deta cut < 2.0)",
0228                                          400,
0229                                          0.,
0230                                          2000.,
0231                                          "M_{jj} WideJets [GeV]",
0232                                          "Number of events");
0233 
0234   m_MjjWide_num_detaL2 = bookH1withSumw2(iBooker,
0235                                          "h1_MjjWide_num_detaL2",
0236                                          "HLT Efficiency Studies (deta cut < 2.0)",
0237                                          400,
0238                                          0.,
0239                                          2000.,
0240                                          "M_{jj} WideJets [GeV]",
0241                                          "Number of events");
0242 
0243   m_MjjWide_den = bookH1withSumw2(iBooker,
0244                                   "h1_MjjWide_den",
0245                                   "HLT Efficiency Studies (default deta cut)",
0246                                   400,
0247                                   0.,
0248                                   2000.,
0249                                   "M_{jj} WideJets [GeV]",
0250                                   "Number of events");
0251 
0252   m_MjjWide_num = bookH1withSumw2(iBooker,
0253                                   "h1_MjjWide_num",
0254                                   "HLT Efficiency Studies (default deta cut)",
0255                                   400,
0256                                   0.,
0257                                   2000.,
0258                                   "M_{jj} WideJets [GeV]",
0259                                   "Number of events");
0260 
0261   m_DetajjWide_finalSel = bookH1withSumw2(iBooker,
0262                                           "h1_DetajjWide_finalSel",
0263                                           "#Delta#eta_{jj} WideJets (final selection)",
0264                                           100,
0265                                           0.,
0266                                           5.,
0267                                           "#Delta#eta_{jj} WideJets",
0268                                           "Number of events");
0269 
0270   m_DetajjWide = bookH1withSumw2(iBooker,
0271                                  "h1_DetajjWide",
0272                                  "#Delta#eta_{jj} WideJets (final selection except #Delta#eta cut)",
0273                                  100,
0274                                  0.,
0275                                  5.,
0276                                  "#Delta#eta_{jj} WideJets",
0277                                  "Number of events");
0278 
0279   m_DphijjWide_finalSel = bookH1withSumw2(iBooker,
0280                                           "h1_DphijjWide_finalSel",
0281                                           "#Delta#phi_{jj} WideJets (final selection)",
0282                                           100,
0283                                           0.,
0284                                           TMath::Pi() + 0.0001,
0285                                           "#Delta#phi_{jj} WideJets [rad.]",
0286                                           "Number of events");
0287 
0288   m_selJets_pt = bookH1withSumw2(
0289       iBooker, "h1_selJets_pt", "Selected CaloJets", 500, 0., 5000., "Jet Pt [GeV]", "Number of events");
0290 
0291   m_selJets_eta =
0292       bookH1withSumw2(iBooker, "h1_selJets_eta", "Selected CaloJets", 100, -5., 5., "#eta", "Number of events");
0293 
0294   m_selJets_phi = bookH1withSumw2(
0295       iBooker, "h1_selJets_phi", "Selected CaloJets", 100, -TMath::Pi(), TMath::Pi(), "#phi (rad.)", "Number of events");
0296 
0297   m_selJets_hadEnergyFraction = bookH1withSumw2(iBooker,
0298                                                 "h1_selJets_hadEnergyFraction",
0299                                                 "Selected CaloJets",
0300                                                 110,
0301                                                 0.,
0302                                                 1.1,
0303                                                 "HAD Energy Fraction",
0304                                                 "Number of events");
0305 
0306   m_selJets_emEnergyFraction = bookH1withSumw2(iBooker,
0307                                                "h1_selJets_emEnergyFraction",
0308                                                "Selected CaloJets",
0309                                                110,
0310                                                0.,
0311                                                1.1,
0312                                                "EM Energy Fraction",
0313                                                "Number of events");
0314 
0315   m_selJets_towersArea = bookH1withSumw2(
0316       iBooker, "h1_selJets_towersArea", "Selected CaloJets", 200, 0., 2., "towers area", "Number of events");
0317 
0318   m_metDiff = bookH1withSumw2(
0319       iBooker, "h1_metDiff", "Met - MetCleaned", 500, -1000., 1000., "met - metcleaned [GeV]", "Number of events");
0320 
0321   m_metCases = bookH1withSumw2(iBooker, "h1_metCases", "Met cases", 3, 0., 3., "case", "Number of events");
0322   m_metCases->getTH1()->GetXaxis()->SetBinLabel(1, "met , metclean");
0323   m_metCases->getTH1()->GetXaxis()->SetBinLabel(2, "met , !metclean");
0324   m_metCases->getTH1()->GetXaxis()->SetBinLabel(3, "!met , metclean");
0325 
0326   m_metCaseNoMetClean = bookH1withSumw2(
0327       iBooker, "h1_metCaseNoMetClean", "Met - MetCleaned", 1000, 0., 2000., "MET [GeV]", "Number of events");
0328 
0329   m_HT_inclusive =
0330       bookH1withSumw2(iBooker, "h1_HT_inclusive", "HT (inclusive)", 150, 0., 15000., "HT [GeV]", "Number of events");
0331 
0332   m_HT_finalSel = bookH1withSumw2(
0333       iBooker, "h1_HT_finalSel", "HT (final selection)", 150, 0., 15000., "HT [GeV]", "Number of events");
0334 
0335   // 2D histograms
0336   m_DetajjVsMjjWide = bookH2withSumw2(iBooker,
0337                                       "h2_DetajjVsMjjWide",
0338                                       "#Delta#eta_{jj} vs M_{jj} WideJets",
0339                                       8000,
0340                                       0.,
0341                                       8000.,
0342                                       100,
0343                                       0.,
0344                                       5.,
0345                                       "M_{jj} WideJets [GeV]",
0346                                       "#Delta#eta_{jj} WideJets");
0347 
0348   m_DetajjVsMjjWide_rebin = bookH2withSumw2(iBooker,
0349                                             "h2_DetajjVsMjjWide_rebin",
0350                                             "#Delta#eta_{jj} vs M_{jj} WideJets",
0351                                             400,
0352                                             0.,
0353                                             8000.,
0354                                             50,
0355                                             0.,
0356                                             5.,
0357                                             "M_{jj} WideJets [GeV]",
0358                                             "#Delta#eta_{jj} WideJets");
0359 
0360   m_metVSmetclean = bookH2withSumw2(
0361       iBooker, "h2_metVSmetclean", "MET clean vs MET", 100, 0., 2000., 100, 0., 2000., "MET [GeV]", "MET clean [GeV]");
0362 }
0363 
0364 // Usual analyze method
0365 void DiJetVarAnalyzer::analyze(const edm::Event &iEvent, const edm::EventSetup &c) {
0366   using namespace std;
0367   using namespace edm;
0368   using namespace reco;
0369 
0370   // ## Get jet collection
0371   edm::Handle<reco::CaloJetCollection> calojets_handle;
0372   iEvent.getByToken(jetCollectionTagToken_, calojets_handle);
0373 
0374   // Loop over all the selected jets ( defined at
0375   // DQM/DataScouting/python/dijetScouting_cff.py )
0376   double thisHT = 0.;
0377   for (reco::CaloJetCollection::const_iterator it = calojets_handle->begin(); it != calojets_handle->end(); ++it) {
0378     // cout << "== jet: " << it->pt() << " " << it->eta() << " " << it->phi() <<
0379     // endl;
0380     m_selJets_pt->Fill(it->pt());
0381     m_selJets_eta->Fill(it->eta());
0382     m_selJets_phi->Fill(it->phi());
0383     m_selJets_hadEnergyFraction->Fill(it->energyFractionHadronic());
0384     m_selJets_emEnergyFraction->Fill(it->emEnergyFraction());
0385     m_selJets_towersArea->Fill(it->towersArea());
0386     thisHT += it->pt();
0387   }
0388 
0389   // HT
0390   m_HT_inclusive->Fill(thisHT);
0391 
0392   // ## Get widejets
0393   edm::Handle<vector<math::PtEtaPhiMLorentzVector>> widejets_handle;
0394   iEvent.getByToken(widejetsCollectionTagToken_, widejets_handle);
0395 
0396   TLorentzVector wj1;
0397   TLorentzVector wj2;
0398   TLorentzVector wdijet;
0399 
0400   double MJJWide = -1;
0401   double DeltaEtaJJWide = -1;
0402   double DeltaPhiJJWide = -1;
0403 
0404   if (widejets_handle->size() >= 2) {
0405     wj1.SetPtEtaPhiM(widejets_handle->at(0).pt(),
0406                      widejets_handle->at(0).eta(),
0407                      widejets_handle->at(0).phi(),
0408                      widejets_handle->at(0).mass());
0409     wj2.SetPtEtaPhiM(widejets_handle->at(1).pt(),
0410                      widejets_handle->at(1).eta(),
0411                      widejets_handle->at(1).phi(),
0412                      widejets_handle->at(1).mass());
0413     wdijet = wj1 + wj2;
0414     MJJWide = wdijet.M();
0415     DeltaEtaJJWide = fabs(wj1.Eta() - wj2.Eta());
0416     DeltaPhiJJWide = fabs(wj1.DeltaPhi(wj2));
0417   }
0418 
0419   // ## Get met collection
0420   // met
0421   edm::Handle<reco::CaloMETCollection> calomet_handle;
0422   iEvent.getByToken(metCollectionTagToken_, calomet_handle);
0423   // met cleaned
0424   edm::Handle<reco::CaloMETCollection> calometClean_handle;
0425   iEvent.getByToken(metCleanCollectionTagToken_, calometClean_handle);
0426 
0427   if (calomet_handle.isValid() && calometClean_handle.isValid()) {
0428     //       std::cout << "---" << std::endl;
0429     //       std::cout << "== calomet: " << (calomet_handle->front()).pt() << "
0430     //       " << (calomet_handle->front()).phi() << std::endl; std::cout << "==
0431     //       calometClean: " << (calometClean_handle->front()).pt() << " " <<
0432     //       (calometClean_handle->front()).phi() << std::endl; std::cout << "==
0433     //       calomet - calometClean: " << (calomet_handle->front()).pt() -
0434     //       (calometClean_handle->front()).pt() << std::endl; std::cout <<
0435     //       "---" << std::endl;
0436     m_metCases->Fill(0);
0437     m_metDiff->Fill((calomet_handle->front()).pt() - (calometClean_handle->front()).pt());
0438     m_metVSmetclean->Fill((calomet_handle->front()).pt(), (calometClean_handle->front()).pt());
0439   } else if (calomet_handle.isValid() && !calometClean_handle.isValid()) {
0440     m_metCases->Fill(1);
0441     m_metCaseNoMetClean->Fill((calomet_handle->front()).pt());
0442   } else if (!calomet_handle.isValid() && calometClean_handle.isValid()) {
0443     m_metCases->Fill(2);
0444   }
0445 
0446   // ## Event Selection
0447   bool pass_nocut = false;
0448   bool pass_twowidejets = false;
0449   bool pass_etaptcuts = false;
0450   bool pass_deta = false;
0451   bool pass_JetIDtwojets = true;
0452   bool pass_dphi = false;
0453   bool pass_metFilter = true;
0454   //--
0455   bool pass_deta_L4 = false;
0456   bool pass_deta_L3 = false;
0457   bool pass_deta_L2 = false;
0458 
0459   bool pass_fullsel_NOdeta = false;
0460   bool pass_fullsel_detaL4 = false;
0461   bool pass_fullsel_detaL3 = false;
0462   bool pass_fullsel_detaL2 = false;
0463   bool pass_fullsel = false;
0464 
0465   // No cut
0466   pass_nocut = true;
0467 
0468   // Two wide jets
0469   if (widejets_handle->size() >= numwidejets_) {
0470     // Two wide jets
0471     pass_twowidejets = true;
0472 
0473     // Eta/pt cuts
0474     if (fabs(wj1.Eta()) < etawidejets_ && wj1.Pt() > ptwidejets_ && fabs(wj2.Eta()) < etawidejets_ &&
0475         wj2.Pt() > ptwidejets_) {
0476       pass_etaptcuts = true;
0477     }
0478 
0479     // Deta cut
0480     if (DeltaEtaJJWide < detawidejets_)
0481       pass_deta = true;
0482 
0483     // Dphi cut
0484     if (DeltaPhiJJWide > dphiwidejets_)
0485       pass_dphi = true;
0486 
0487     // Other Deta cuts
0488     if (DeltaEtaJJWide < 4.0)
0489       pass_deta_L4 = true;
0490 
0491     if (DeltaEtaJJWide < 3.0)
0492       pass_deta_L3 = true;
0493 
0494     if (DeltaEtaJJWide < 2.0)
0495       pass_deta_L2 = true;
0496   }
0497   // Jet id two leading jets
0498   if (calojets_handle->size() >= numwidejets_) {
0499     //   first jet
0500     reco::CaloJetCollection::const_iterator thisJet = calojets_handle->begin();
0501     // cout << "== thisJet1: " << thisJet->pt() << " " << thisJet->eta() << " "
0502     // << thisJet->phi() << endl;
0503     if (thisJet->energyFractionHadronic() > maxHADfraction_ || thisJet->emEnergyFraction() > maxEMfraction_)
0504       pass_JetIDtwojets = false;
0505 
0506     //   second jet
0507     ++thisJet;
0508     // cout << "== thisJet2: " << thisJet->pt() << " " << thisJet->eta() << " "
0509     // << thisJet->phi() << endl;
0510     if (thisJet->energyFractionHadronic() > maxHADfraction_ || thisJet->emEnergyFraction() > maxEMfraction_)
0511       pass_JetIDtwojets = false;
0512   }
0513   // Met filter
0514   if (calomet_handle.isValid() && calometClean_handle.isValid()) {
0515     if (fabs((calomet_handle->front()).pt() - (calometClean_handle->front()).pt()) > 0.1)
0516       pass_metFilter = false;
0517   }
0518   // Full selection (no deta cut)
0519   if (pass_nocut && pass_twowidejets && pass_etaptcuts && pass_JetIDtwojets && pass_dphi && pass_metFilter)
0520     pass_fullsel_NOdeta = true;
0521   // Full selection (various deta cuts)
0522   if (pass_nocut && pass_twowidejets && pass_etaptcuts && pass_JetIDtwojets && pass_dphi && pass_metFilter &&
0523       pass_deta_L4)
0524     pass_fullsel_detaL4 = true;
0525   if (pass_nocut && pass_twowidejets && pass_etaptcuts && pass_JetIDtwojets && pass_dphi && pass_metFilter &&
0526       pass_deta_L3)
0527     pass_fullsel_detaL3 = true;
0528   if (pass_nocut && pass_twowidejets && pass_etaptcuts && pass_JetIDtwojets && pass_dphi && pass_metFilter &&
0529       pass_deta_L2)
0530     pass_fullsel_detaL2 = true;
0531   // Full selection (default deta cut)
0532   if (pass_nocut && pass_twowidejets && pass_etaptcuts && pass_deta && pass_JetIDtwojets && pass_dphi && pass_metFilter)
0533     pass_fullsel = true;
0534 
0535   // ## Fill Histograms
0536   // Cut-flow plot
0537   if (pass_nocut)
0538     m_cutFlow->Fill(0);
0539   if (pass_nocut && pass_twowidejets)
0540     m_cutFlow->Fill(1);
0541   if (pass_nocut && pass_twowidejets && pass_etaptcuts)
0542     m_cutFlow->Fill(2);
0543   if (pass_nocut && pass_twowidejets && pass_etaptcuts && pass_deta)
0544     m_cutFlow->Fill(3);
0545   if (pass_nocut && pass_twowidejets && pass_etaptcuts && pass_deta && pass_JetIDtwojets)
0546     m_cutFlow->Fill(4);
0547   if (pass_nocut && pass_twowidejets && pass_etaptcuts && pass_deta && pass_JetIDtwojets && pass_dphi)
0548     m_cutFlow->Fill(5);
0549   if (pass_fullsel)
0550     m_cutFlow->Fill(6);
0551 
0552   // After full selection
0553   if (pass_fullsel) {
0554     // 1D histograms
0555     m_MjjWide_finalSel->Fill(MJJWide);
0556     m_MjjWide_finalSel_varbin->Fill(MJJWide);
0557     m_DetajjWide_finalSel->Fill(DeltaEtaJJWide);
0558     m_DphijjWide_finalSel->Fill(DeltaPhiJJWide);
0559     m_HT_finalSel->Fill(thisHT);
0560   }
0561   // After full selection (without "noise" filters)
0562   if (pass_nocut && pass_twowidejets && pass_etaptcuts && pass_deta) {
0563     m_MjjWide_finalSel_WithoutNoiseFilter->Fill(MJJWide);
0564     m_MjjWide_finalSel_WithoutNoiseFilter_varbin->Fill(MJJWide);
0565   }
0566   // After full selection (except DeltaEta cut)
0567   if (pass_fullsel_NOdeta) {
0568     // 1D histograms
0569     m_DetajjWide->Fill(DeltaEtaJJWide);
0570     if (DeltaEtaJJWide >= 0.0 && DeltaEtaJJWide < 0.5)
0571       m_MjjWide_deta_0p0_0p5->Fill(MJJWide);
0572     if (DeltaEtaJJWide >= 0.5 && DeltaEtaJJWide < 1.0)
0573       m_MjjWide_deta_0p5_1p0->Fill(MJJWide);
0574     if (DeltaEtaJJWide >= 1.0 && DeltaEtaJJWide < 1.5)
0575       m_MjjWide_deta_1p0_1p5->Fill(MJJWide);
0576     if (DeltaEtaJJWide >= 1.5 && DeltaEtaJJWide < 2.0)
0577       m_MjjWide_deta_1p5_2p0->Fill(MJJWide);
0578     if (DeltaEtaJJWide >= 2.0 && DeltaEtaJJWide < 2.5)
0579       m_MjjWide_deta_2p0_2p5->Fill(MJJWide);
0580     if (DeltaEtaJJWide >= 2.5 && DeltaEtaJJWide < 3.0)
0581       m_MjjWide_deta_2p5_3p0->Fill(MJJWide);
0582     if (DeltaEtaJJWide >= 3.0)
0583       m_MjjWide_deta_3p0_inf->Fill(MJJWide);
0584 
0585     // 2D histograms
0586     m_DetajjVsMjjWide->Fill(MJJWide, DeltaEtaJJWide);
0587     m_DetajjVsMjjWide_rebin->Fill(MJJWide, DeltaEtaJJWide);
0588   }
0589 
0590   // ## Get Trigger Info
0591 
0592   // HLT paths for DataScouting
0593   //  DST_HT250_v1
0594   //  DST_L1HTT_Or_L1MultiJet_v1
0595   //  DST_Mu5_HT250_v1
0596   //  DST_Ele8_CaloIdL_CaloIsoVL_TrkIdVL_TrkIsoVL_HT250_v1
0597   int HLTpathMain_fired = -1;
0598   int HLTpathMonitor_fired = -1;
0599   if (HLTpathMain_ and HLTpathMonitor_ and triggerConfiguration_.setEvent(iEvent, c)) {
0600     // invalid HLT configuration, skip the processing
0601 
0602     // if the L1 or HLT configurations have changed, (re)initialize the filters
0603     // (including during the first event)
0604     if (triggerConfiguration_.configurationUpdated()) {
0605       HLTpathMain_->init(triggerConfiguration_);
0606       HLTpathMonitor_->init(triggerConfiguration_);
0607 
0608       // log the expanded configuration
0609       // std::cout << "HLT selector configurations updated" << std::endl;
0610       // std::cout << "HLTpathMain:    " << *HLTpathMain_    << std::endl;
0611       // std::cout << "HLTpathMonitor: " << *HLTpathMonitor_ << std::endl;
0612     }
0613 
0614     HLTpathMain_fired = (*HLTpathMain_)(triggerConfiguration_);
0615     HLTpathMonitor_fired = (*HLTpathMonitor_)(triggerConfiguration_);
0616 
0617     // The OR of the two should always be "1"
0618     // std::cout << *HLTpathMain_ << ": " << HLTpathMain_fired << " -- " <<
0619     // *HLTpathMonitor_ << ": " << HLTpathMonitor_fired << std::endl;
0620   }
0621 
0622   // ## Trigger Efficiency Curves
0623 
0624   // denominator - full sel NO deta cut
0625   if (pass_fullsel_NOdeta && HLTpathMonitor_fired == 1) {
0626     m_MjjWide_den_NOdeta->Fill(MJJWide);
0627 
0628     // numerator
0629     if (HLTpathMain_fired == 1) {
0630       m_MjjWide_num_NOdeta->Fill(MJJWide);
0631     }
0632   }
0633 
0634   // denominator - full sel deta < 4.0
0635   if (pass_fullsel_detaL4 && HLTpathMonitor_fired == 1) {
0636     m_MjjWide_den_detaL4->Fill(MJJWide);
0637 
0638     // numerator
0639     if (HLTpathMain_fired == 1) {
0640       m_MjjWide_num_detaL4->Fill(MJJWide);
0641     }
0642   }
0643 
0644   // denominator - full sel deta < 3.0
0645   if (pass_fullsel_detaL3 && HLTpathMonitor_fired == 1) {
0646     m_MjjWide_den_detaL3->Fill(MJJWide);
0647 
0648     // numerator
0649     if (HLTpathMain_fired == 1) {
0650       m_MjjWide_num_detaL3->Fill(MJJWide);
0651     }
0652   }
0653 
0654   // denominator - full sel deta < 2.0
0655   if (pass_fullsel_detaL2 && HLTpathMonitor_fired == 1) {
0656     m_MjjWide_den_detaL2->Fill(MJJWide);
0657 
0658     // numerator
0659     if (HLTpathMain_fired == 1) {
0660       m_MjjWide_num_detaL2->Fill(MJJWide);
0661     }
0662   }
0663 
0664   // denominator - full sel default deta cut (typically 1.3)
0665   if (pass_fullsel && HLTpathMonitor_fired == 1) {
0666     m_MjjWide_den->Fill(MJJWide);
0667 
0668     // numerator
0669     if (HLTpathMain_fired == 1) {
0670       m_MjjWide_num->Fill(MJJWide);
0671     }
0672   }
0673 }