Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:18:57

0001 #include "DataFormats/HLTReco/interface/TriggerObject.h"
0002 #include "FWCore/Common/interface/TriggerNames.h"
0003 #include "FWCore/Framework/interface/Frameworkfwd.h"
0004 #include "FWCore/Framework/interface/MakerMacros.h"
0005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0006 #include "HLTriggerOffline/SUSYBSM/interface/SUSY_HLT_alphaT.h"
0007 //#include "HLTriggerOffline/SUSYBSM/interface/AlphaT.h"
0008 
0009 typedef ROOT::Math::LorentzVector<ROOT::Math::PtEtaPhiM4D<double>> LorentzV;
0010 
0011 SUSY_HLT_alphaT::SUSY_HLT_alphaT(const edm::ParameterSet &ps) {
0012   edm::LogInfo("SUSY_HLT_alphaT") << "Constructor SUSY_HLT_alphaT::SUSY_HLT_alphaT " << std::endl;
0013   // Get parameters from configuration file
0014   theTrigSummary_ = consumes<trigger::TriggerEvent>(ps.getParameter<edm::InputTag>("trigSummary"));
0015   thePfJetCollection_ = consumes<reco::PFJetCollection>(ps.getParameter<edm::InputTag>("pfJetCollection"));
0016   // theCaloJetCollection_ =
0017   // consumes<reco::CaloJetCollection>(ps.getParameter<edm::InputTag>("caloJetCollection"));
0018   triggerResults_ = consumes<edm::TriggerResults>(ps.getParameter<edm::InputTag>("TriggerResults"));
0019   HLTProcess_ = ps.getParameter<std::string>("HLTProcess");
0020   triggerPath_ = ps.getParameter<std::string>("TriggerPath");
0021   triggerPathAuxiliaryForHadronic_ = ps.getParameter<std::string>("TriggerPathAuxiliaryForHadronic");
0022   triggerFilter_ = ps.getParameter<edm::InputTag>("TriggerFilter");
0023   triggerPreFilter_ = ps.getParameter<edm::InputTag>("TriggerPreFilter");
0024   ptThrJet_ = ps.getUntrackedParameter<double>("PtThrJet");
0025   etaThrJet_ = ps.getUntrackedParameter<double>("EtaThrJet");
0026   // caloAlphaTThrTurnon_ =
0027   // ps.getUntrackedParameter<double>("caloAlphaTThrTurnon"); caloHtThrTurnon_ =
0028   // ps.getUntrackedParameter<double>("caloHtThrTurnon");
0029   pfAlphaTThrTurnon_ = ps.getUntrackedParameter<double>("pfAlphaTThrTurnon");
0030   pfHtThrTurnon_ = ps.getUntrackedParameter<double>("pfHtThrTurnon");
0031 }
0032 
0033 SUSY_HLT_alphaT::~SUSY_HLT_alphaT() {
0034   edm::LogInfo("SUSY_HLT_alphaT") << "Destructor SUSY_HLT_alphaT::~SUSY_HLT_alphaT " << std::endl;
0035 }
0036 
0037 void SUSY_HLT_alphaT::dqmBeginRun(edm::Run const &run, edm::EventSetup const &e) {
0038   bool changed;
0039 
0040   if (!fHltConfig.init(run, e, HLTProcess_, changed)) {
0041     edm::LogError("SUSY_HLT_alphaT") << "Initialization of HLTConfigProvider failed!!";
0042     return;
0043   }
0044 
0045   bool pathFound = false;
0046   const std::vector<std::string> allTrigNames = fHltConfig.triggerNames();
0047   for (size_t j = 0; j < allTrigNames.size(); ++j) {
0048     if (allTrigNames[j].find(triggerPath_) != std::string::npos) {
0049       pathFound = true;
0050     }
0051   }
0052 
0053   if (!pathFound) {
0054     LogDebug("SUSY_HLT_alphaT") << "Path not found"
0055                                 << "\n";
0056     return;
0057   }
0058   // std::vector<std::string> filtertags = fHltConfig.moduleLabels( triggerPath_
0059   // ); triggerFilter_ =
0060   // edm::InputTag(filtertags[filtertags.size()-1],"",fHltConfig.processName());
0061   // triggerFilter_ = edm::InputTag("hltPFMET120Mu5L3PreFiltered", "",
0062   // fHltConfig.processName());
0063 
0064   edm::LogInfo("SUSY_HLT_alphaT") << "SUSY_HLT_alphaT::beginRun" << std::endl;
0065 }
0066 
0067 void SUSY_HLT_alphaT::bookHistograms(DQMStore::IBooker &ibooker_, edm::Run const &, edm::EventSetup const &) {
0068   edm::LogInfo("SUSY_HLT_alphaT") << "SUSY_HLT_alphaT::bookHistograms" << std::endl;
0069   // book at beginRun
0070   bookHistos(ibooker_);
0071 }
0072 
0073 void SUSY_HLT_alphaT::analyze(edm::Event const &e, edm::EventSetup const &eSetup) {
0074   edm::LogInfo("SUSY_HLT_alphaT") << "SUSY_HLT_alphaT::analyze" << std::endl;
0075 
0076   //-------------------------------
0077   //--- Trigger
0078   //-------------------------------
0079   edm::Handle<edm::TriggerResults> hltresults;
0080   e.getByToken(triggerResults_, hltresults);
0081   if (!hltresults.isValid()) {
0082     edm::LogWarning("SUSY_HLT_alphaT") << "invalid collection: TriggerResults"
0083                                        << "\n";
0084     return;
0085   }
0086   edm::Handle<trigger::TriggerEvent> triggerSummary;
0087   e.getByToken(theTrigSummary_, triggerSummary);
0088   if (!triggerSummary.isValid()) {
0089     edm::LogWarning("SUSY_HLT_alphaT") << "invalid collection: TriggerSummary"
0090                                        << "\n";
0091     return;
0092   }
0093 
0094   //-------------------------------
0095   //--- Jets
0096   //-------------------------------
0097   edm::Handle<reco::PFJetCollection> pfJetCollection;
0098   e.getByToken(thePfJetCollection_, pfJetCollection);
0099   if (!pfJetCollection.isValid()) {
0100     edm::LogWarning("SUSY_HLT_alphaT") << "invalid collection: PFJets"
0101                                        << "\n";
0102     return;
0103   }
0104   // edm::Handle<reco::CaloJetCollection> caloJetCollection;
0105   // e.getByToken (theCaloJetCollection_,caloJetCollection);
0106   // if ( !caloJetCollection.isValid() ){
0107   //     edm::LogWarning ("SUSY_HLT_alphaT") << "invalid collection: CaloJets"
0108   //     << "\n"; return;
0109   // }
0110 
0111   // get online objects
0112   // For now just get the jets and recalculate ht and alphaT
0113   size_t filterIndex = triggerSummary->filterIndex(triggerFilter_);
0114   // size_t preFilterIndex = triggerSummary->filterIndex( triggerPreFilter_ );
0115   trigger::TriggerObjectCollection triggerObjects = triggerSummary->getObjects();
0116 
0117   // Get the PF objects from the filter
0118   double hltPfHt = 0.;
0119   std::vector<LorentzV> hltPfJets;
0120   if (!(filterIndex >= triggerSummary->sizeFilters())) {
0121     const trigger::Keys &keys = triggerSummary->filterKeys(filterIndex);
0122 
0123     for (size_t j = 0; j < keys.size(); ++j) {
0124       trigger::TriggerObject foundObject = triggerObjects[keys[j]];
0125 
0126       //  if(foundObject.id() == 85){ //It's a jet
0127       if (foundObject.pt() > ptThrJet_ && fabs(foundObject.eta()) < etaThrJet_) {
0128         hltPfHt += foundObject.pt();
0129         LorentzV JetLVec(foundObject.pt(), foundObject.eta(), foundObject.phi(), foundObject.mass());
0130         hltPfJets.push_back(JetLVec);
0131       }
0132       //   }
0133     }
0134   }
0135 
0136   // //Get the Calo objects from the prefilter
0137   // double hltCaloHt=0.;
0138   // std::vector<LorentzV> hltCaloJets;
0139   // if( !(preFilterIndex >= triggerSummary->sizeFilters()) ){
0140   //     const trigger::Keys& keys = triggerSummary->filterKeys( preFilterIndex
0141   //     );
0142 
0143   //     for( size_t j = 0; j < keys.size(); ++j ){
0144   //         trigger::TriggerObject foundObject = triggerObjects[keys[j]];
0145 
0146   //         //  if(foundObject.id() == 85){ //It's a jet
0147   //         if(foundObject.pt()>ptThrJet_ && fabs(foundObject.eta()) <
0148   //         etaThrJet_){
0149   //             hltCaloHt += foundObject.pt();
0150   //             LorentzV
0151   //             JetLVec(foundObject.pt(),foundObject.eta(),foundObject.phi(),foundObject.mass());
0152   //             hltCaloJets.push_back(JetLVec);
0153   //         }
0154   //         //   }
0155   //     }
0156   // }
0157 
0158   // Fill the alphaT and HT histograms
0159   if (!hltPfJets.empty()) {
0160     double hltPfAlphaT = AlphaT(hltPfJets, true).value();
0161     h_triggerPfAlphaT->Fill(hltPfAlphaT);
0162     h_triggerPfHt->Fill(hltPfHt);
0163     h_triggerPfAlphaT_triggerPfHt->Fill(hltPfHt, hltPfAlphaT);
0164   }
0165 
0166   // if(hltCaloJets.size()>0){
0167   //     double hltCaloAlphaT = AlphaT(hltCaloJets,true).value();
0168   //     h_triggerCaloAlphaT->Fill(hltCaloAlphaT);
0169   //     h_triggerCaloHt->Fill(hltCaloHt);
0170   //     h_triggerCaloAlphaT_triggerCaloHt->Fill(hltCaloHt, hltCaloAlphaT);
0171   // }
0172 
0173   bool hasFired = false;
0174   bool hasFiredAuxiliaryForHadronicLeg = false;
0175   const edm::TriggerNames &trigNames = e.triggerNames(*hltresults);
0176   unsigned int numTriggers = trigNames.size();
0177   for (unsigned int hltIndex = 0; hltIndex < numTriggers; ++hltIndex) {
0178     if (trigNames.triggerName(hltIndex).find(triggerPath_) != std::string::npos && hltresults->wasrun(hltIndex) &&
0179         hltresults->accept(hltIndex))
0180       hasFired = true;
0181     if (trigNames.triggerName(hltIndex).find(triggerPathAuxiliaryForHadronic_) != std::string::npos &&
0182         hltresults->wasrun(hltIndex) && hltresults->accept(hltIndex))
0183       hasFiredAuxiliaryForHadronicLeg = true;
0184   }
0185 
0186   if (hasFiredAuxiliaryForHadronicLeg) {
0187     float pfHT = 0.0;
0188     std::vector<LorentzV> pfJets;
0189     for (reco::PFJetCollection::const_iterator i_pfjet = pfJetCollection->begin(); i_pfjet != pfJetCollection->end();
0190          ++i_pfjet) {
0191       if (i_pfjet->pt() < ptThrJet_)
0192         continue;
0193       if (fabs(i_pfjet->eta()) > etaThrJet_)
0194         continue;
0195       pfHT += i_pfjet->pt();
0196       LorentzV JetLVec(i_pfjet->pt(), i_pfjet->eta(), i_pfjet->phi(), i_pfjet->mass());
0197       pfJets.push_back(JetLVec);
0198     }
0199 
0200     // //Make the gen Calo HT and AlphaT
0201     // float caloHT = 0.0;
0202     // std::vector<LorentzV> caloJets;
0203     // for (reco::CaloJetCollection::const_iterator i_calojet =
0204     // caloJetCollection->begin(); i_calojet != caloJetCollection->end();
0205     // ++i_calojet){
0206     //   if (i_calojet->pt() < ptThrJet_) continue;
0207     //   if (fabs(i_calojet->eta()) > etaThrJet_) continue;
0208     //   caloHT += i_calojet->pt();
0209     //   LorentzV
0210     //   JetLVec(i_calojet->pt(),i_calojet->eta(),i_calojet->phi(),i_calojet->mass());
0211     //   caloJets.push_back(JetLVec);
0212     // }
0213 
0214     // AlphaT aT = AlphaT(jets);
0215     // double caloAlphaT = AlphaT(caloJets).value();
0216     double pfAlphaT = AlphaT(pfJets).value();
0217 
0218     // Fill the turnons
0219     if (hasFired) {
0220       // if(caloHT>caloHtThrTurnon_) h_caloAlphaTTurnOn_num-> Fill(caloAlphaT);
0221       // if(caloAlphaT>caloAlphaTThrTurnon_) h_caloHtTurnOn_num-> Fill(caloHT);
0222 
0223       if (pfHT > pfHtThrTurnon_)
0224         h_pfAlphaTTurnOn_num->Fill(pfAlphaT);
0225       if (pfAlphaT > pfAlphaTThrTurnon_)
0226         h_pfHtTurnOn_num->Fill(pfHT);
0227     }
0228     // if(caloHT>caloHtThrTurnon_) h_caloAlphaTTurnOn_den-> Fill(caloAlphaT);
0229     // if(caloAlphaT>caloAlphaTThrTurnon_) h_caloHtTurnOn_den-> Fill(caloHT);
0230 
0231     if (pfHT > pfHtThrTurnon_)
0232       h_pfAlphaTTurnOn_den->Fill(pfAlphaT);
0233     if (pfAlphaT > pfAlphaTThrTurnon_)
0234       h_pfHtTurnOn_den->Fill(pfHT);
0235   }
0236 }
0237 
0238 void SUSY_HLT_alphaT::bookHistos(DQMStore::IBooker &ibooker_) {
0239   ibooker_.cd();
0240   ibooker_.setCurrentFolder("HLT/SUSYBSM/" + triggerPath_);
0241 
0242   // offline quantities
0243 
0244   // online quantities
0245   // h_triggerCaloHt = ibooker_.book1D("triggerCaloHt", "Trigger Calo Ht; HT
0246   // (GeV)", 60, 0.0, 1500.0); h_triggerCaloAlphaT =
0247   // ibooker_.book1D("triggerCaloAlphaT", "Trigger Calo AlphaT; AlphaT", 80,
0248   // 0., 1.0); h_triggerCaloAlphaT_triggerCaloHt =
0249   // ibooker_.book2D("triggerCaloAlphaT_triggerCaloHt","Trigger Calo HT vs
0250   // Trigger Calo AlphaT; HT (GeV); AlphaT", 60,0.0,1500.,80,0.,1.0);
0251   h_triggerPfHt = ibooker_.book1D("triggerPfHt", "Trigger PF Ht; HT (GeV)", 60, 0.0, 1500.0);
0252   h_triggerPfAlphaT = ibooker_.book1D("triggerPfAlphaT", "Trigger PF AlphaT; AlphaT", 80, 0., 1.0);
0253   h_triggerPfAlphaT_triggerPfHt = ibooker_.book2D("triggerPfAlphaT_triggerPfHt",
0254                                                   "Trigger PF HT vs Trigger PF AlphaT; HT (GeV); AlphaT",
0255                                                   60,
0256                                                   0.0,
0257                                                   1500.,
0258                                                   80,
0259                                                   0.,
0260                                                   1.0);
0261 
0262   // num and den hists to be divided in harvesting step to make turn on curves
0263   // h_caloAlphaTTurnOn_num = ibooker_.book1D("caloAlphaTTurnOn_num", "Calo
0264   // AlphaT Turn On Numerator; AlphaT", 40, 0.0, 1.0 ); h_caloAlphaTTurnOn_den =
0265   // ibooker_.book1D("caloAlphaTTurnOn_den", "Calo AlphaT Turn OnDenominator;
0266   // AlphaT", 40, 0.0, 1.0 ); h_caloHtTurnOn_num =
0267   // ibooker_.book1D("caloHtTurnOn_num", "Calo HT Turn On Numerator; HT (GeV)",
0268   // 30, 0.0, 1500.0 ); h_caloHtTurnOn_den = ibooker_.book1D("caloHtTurnOn_den",
0269   // "Calo HT Turn On Denominator; HT (GeV)", 30, 0.0, 1500.0 );
0270 
0271   h_pfAlphaTTurnOn_num = ibooker_.book1D("pfAlphaTTurnOn_num", "PF AlphaT Turn On Numerator; AlphaT", 40, 0.0, 1.0);
0272   h_pfAlphaTTurnOn_den = ibooker_.book1D("pfAlphaTTurnOn_den", "PF AlphaT Turn OnDenominator; AlphaT", 40, 0.0, 1.0);
0273   h_pfHtTurnOn_num = ibooker_.book1D("pfHtTurnOn_num", "PF HT Turn On Numerator; HT (GeV)", 30, 0.0, 1500.0);
0274   h_pfHtTurnOn_den = ibooker_.book1D("pfHtTurnOn_den", "PF HT Turn On Denominator; HT (GeV)", 30, 0.0, 1500.0);
0275 
0276   ibooker_.cd();
0277 }
0278 
0279 // define this as a plug-in
0280 DEFINE_FWK_MODULE(SUSY_HLT_alphaT);