Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:10:07

0001 #include "DataFormats/HLTReco/interface/TriggerObject.h"
0002 #include "DataFormats/Math/interface/LorentzVector.h"
0003 #include "FWCore/Common/interface/TriggerNames.h"
0004 #include "FWCore/Framework/interface/Frameworkfwd.h"
0005 #include "FWCore/Framework/interface/MakerMacros.h"
0006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0007 #include "HLTriggerOffline/SUSYBSM/interface/SUSY_HLT_DoubleEle_Hadronic.h"
0008 
0009 SUSY_HLT_DoubleEle_Hadronic::SUSY_HLT_DoubleEle_Hadronic(const edm::ParameterSet &ps) {
0010   edm::LogInfo("SUSY_HLT_DoubleEle_Hadronic")
0011       << "Constructor SUSY_HLT_DoubleEle_Hadronic::SUSY_HLT_DoubleEle_Hadronic " << std::endl;
0012   // Get parameters from configuration file

0013   theTrigSummary_ = consumes<trigger::TriggerEvent>(ps.getParameter<edm::InputTag>("trigSummary"));
0014   theElectronCollection_ = consumes<reco::GsfElectronCollection>(ps.getParameter<edm::InputTag>("ElectronCollection"));
0015   thePfJetCollection_ = consumes<reco::PFJetCollection>(ps.getParameter<edm::InputTag>("pfJetCollection"));
0016   theCaloJetCollection_ = consumes<reco::CaloJetCollection>(ps.getParameter<edm::InputTag>("caloJetCollection"));
0017   triggerResults_ = consumes<edm::TriggerResults>(ps.getParameter<edm::InputTag>("TriggerResults"));
0018   HLTProcess_ = ps.getParameter<std::string>("HLTProcess");
0019   triggerPath_ = ps.getParameter<std::string>("TriggerPath");
0020   triggerPathAuxiliaryForElectron_ = ps.getParameter<std::string>("TriggerPathAuxiliaryForElectron");
0021   triggerPathAuxiliaryForHadronic_ = ps.getParameter<std::string>("TriggerPathAuxiliaryForHadronic");
0022   triggerFilter_ = ps.getParameter<edm::InputTag>("TriggerFilter");
0023   ptThrJet_ = ps.getUntrackedParameter<double>("PtThrJet");
0024   etaThrJet_ = ps.getUntrackedParameter<double>("EtaThrJet");
0025 }
0026 
0027 SUSY_HLT_DoubleEle_Hadronic::~SUSY_HLT_DoubleEle_Hadronic() {
0028   edm::LogInfo("SUSY_HLT_DoubleEle_Hadronic")
0029       << "Destructor SUSY_HLT_DoubleEle_Hadronic::~SUSY_HLT_DoubleEle_Hadronic " << std::endl;
0030 }
0031 
0032 void SUSY_HLT_DoubleEle_Hadronic::dqmBeginRun(edm::Run const &run, edm::EventSetup const &e) {
0033   bool changed;
0034 
0035   if (!fHltConfig.init(run, e, HLTProcess_, changed)) {
0036     edm::LogError("SUSY_HLT_DoubleEle_Hadronic") << "Initialization of HLTConfigProvider failed!!";
0037     return;
0038   }
0039 
0040   bool pathFound = false;
0041   const std::vector<std::string> allTrigNames = fHltConfig.triggerNames();
0042   for (size_t j = 0; j < allTrigNames.size(); ++j) {
0043     if (allTrigNames[j].find(triggerPath_) != std::string::npos) {
0044       pathFound = true;
0045     }
0046   }
0047 
0048   if (!pathFound) {
0049     LogDebug("SUSY_HLT_DoubleEle_Hadronic") << "Path not found"
0050                                             << "\n";
0051     return;
0052   }
0053   // std::vector<std::string> filtertags = fHltConfig.moduleLabels( triggerPath_

0054   // ); triggerFilter_ =

0055   // edm::InputTag(filtertags[filtertags.size()-1],"",fHltConfig.processName());

0056   // triggerFilter_ = edm::InputTag("hltPFMET120Mu5L3PreFiltered", "",

0057   // fHltConfig.processName());

0058 
0059   edm::LogInfo("SUSY_HLT_DoubleEle_Hadronic") << "SUSY_HLT_DoubleEle_Hadronic::beginRun" << std::endl;
0060 }
0061 
0062 void SUSY_HLT_DoubleEle_Hadronic::bookHistograms(DQMStore::IBooker &ibooker_,
0063                                                  edm::Run const &,
0064                                                  edm::EventSetup const &) {
0065   edm::LogInfo("SUSY_HLT_DoubleEle_Hadronic") << "SUSY_HLT_DoubleEle_Hadronic::bookHistograms" << std::endl;
0066   // book at beginRun

0067   bookHistos(ibooker_);
0068 }
0069 
0070 void SUSY_HLT_DoubleEle_Hadronic::analyze(edm::Event const &e, edm::EventSetup const &eSetup) {
0071   edm::LogInfo("SUSY_HLT_DoubleEle_Hadronic") << "SUSY_HLT_DoubleEle_Hadronic::analyze" << std::endl;
0072 
0073   //-------------------------------

0074   //--- Jets

0075   //-------------------------------

0076   edm::Handle<reco::PFJetCollection> pfJetCollection;
0077   e.getByToken(thePfJetCollection_, pfJetCollection);
0078   if (!pfJetCollection.isValid()) {
0079     edm::LogError("SUSY_HLT_DoubleEle_Hadronic") << "invalid collection: PFJets"
0080                                                  << "\n";
0081     return;
0082   }
0083   edm::Handle<reco::CaloJetCollection> caloJetCollection;
0084   e.getByToken(theCaloJetCollection_, caloJetCollection);
0085   if (!caloJetCollection.isValid()) {
0086     edm::LogError("SUSY_HLT_DoubleEle_Hadronic") << "invalid collection: CaloJets"
0087                                                  << "\n";
0088     return;
0089   }
0090 
0091   //-------------------------------

0092   //--- Electron

0093   //-------------------------------

0094   edm::Handle<reco::GsfElectronCollection> ElectronCollection;
0095   e.getByToken(theElectronCollection_, ElectronCollection);
0096   if (!ElectronCollection.isValid()) {
0097     edm::LogError("SUSY_HLT_DoubleEle_Hadronic") << "invalid collection: Electrons "
0098                                                  << "\n";
0099     return;
0100   }
0101 
0102   //-------------------------------

0103   //--- Trigger

0104   //-------------------------------

0105   edm::Handle<edm::TriggerResults> hltresults;
0106   e.getByToken(triggerResults_, hltresults);
0107   if (!hltresults.isValid()) {
0108     edm::LogError("SUSY_HLT_DoubleEle_Hadronic") << "invalid collection: TriggerResults"
0109                                                  << "\n";
0110     return;
0111   }
0112   edm::Handle<trigger::TriggerEvent> triggerSummary;
0113   e.getByToken(theTrigSummary_, triggerSummary);
0114   if (!triggerSummary.isValid()) {
0115     edm::LogError("SUSY_HLT_DoubleEle_Hadronic") << "invalid collection: TriggerSummary"
0116                                                  << "\n";
0117     return;
0118   }
0119 
0120   // get online objects

0121   std::vector<float> ptElectron, etaElectron, phiElectron;
0122   size_t filterIndex = triggerSummary->filterIndex(triggerFilter_);
0123   trigger::TriggerObjectCollection triggerObjects = triggerSummary->getObjects();
0124   if (!(filterIndex >= triggerSummary->sizeFilters())) {
0125     const trigger::Keys &keys = triggerSummary->filterKeys(filterIndex);
0126     for (size_t j = 0; j < keys.size(); ++j) {
0127       trigger::TriggerObject foundObject = triggerObjects[keys[j]];
0128       if (fabs(foundObject.id()) == 11) {  // It's an electron

0129 
0130         bool same = false;
0131         for (unsigned int x = 0; x < ptElectron.size(); x++) {
0132           if (fabs(ptElectron[x] - foundObject.pt()) < 0.01 || fabs(etaElectron[x] - foundObject.eta()) < 0.001 ||
0133               fabs(phiElectron[x] - foundObject.phi()) < 0.001)
0134             same = true;
0135         }
0136 
0137         if (!same) {
0138           h_triggerElePt->Fill(foundObject.pt());
0139           h_triggerEleEta->Fill(foundObject.eta());
0140           h_triggerElePhi->Fill(foundObject.phi());
0141           ptElectron.push_back(foundObject.pt());
0142           etaElectron.push_back(foundObject.eta());
0143           phiElectron.push_back(foundObject.phi());
0144         }
0145       }
0146     }
0147     if (ptElectron.size() >= 2) {
0148       math::PtEtaPhiMLorentzVectorD *ele1 =
0149           new math::PtEtaPhiMLorentzVectorD(ptElectron[0], etaElectron[0], phiElectron[0], 0.0005);
0150       math::PtEtaPhiMLorentzVectorD *ele2 =
0151           new math::PtEtaPhiMLorentzVectorD(ptElectron[1], etaElectron[1], phiElectron[1], 0.0005);
0152       (*ele1) += (*ele2);
0153       h_triggerDoubleEleMass->Fill(ele1->M());
0154       delete ele1;
0155       delete ele2;
0156     } else {
0157       h_triggerDoubleEleMass->Fill(-1.);
0158     }
0159   }
0160 
0161   bool hasFired = false;
0162   bool hasFiredAuxiliaryForElectronLeg = false;
0163   bool hasFiredAuxiliaryForHadronicLeg = false;
0164   const edm::TriggerNames &trigNames = e.triggerNames(*hltresults);
0165   unsigned int numTriggers = trigNames.size();
0166   for (unsigned int hltIndex = 0; hltIndex < numTriggers; ++hltIndex) {
0167     if (trigNames.triggerName(hltIndex).find(triggerPath_) != std::string::npos && hltresults->wasrun(hltIndex) &&
0168         hltresults->accept(hltIndex))
0169       hasFired = true;
0170     if (trigNames.triggerName(hltIndex).find(triggerPathAuxiliaryForElectron_) != std::string::npos &&
0171         hltresults->wasrun(hltIndex) && hltresults->accept(hltIndex))
0172       hasFiredAuxiliaryForElectronLeg = true;
0173     if (trigNames.triggerName(hltIndex).find(triggerPathAuxiliaryForHadronic_) != std::string::npos &&
0174         hltresults->wasrun(hltIndex) && hltresults->accept(hltIndex))
0175       hasFiredAuxiliaryForHadronicLeg = true;
0176   }
0177 
0178   if (hasFiredAuxiliaryForElectronLeg || hasFiredAuxiliaryForHadronicLeg) {
0179     // Matching the Electron

0180     int indexOfMatchedElectron[2] = {-1};
0181     int matchedCounter = 0;
0182     int offlineCounter = 0;
0183     for (reco::GsfElectronCollection::const_iterator Electron = ElectronCollection->begin();
0184          (Electron != ElectronCollection->end() && matchedCounter < 2);
0185          ++Electron) {
0186       for (size_t off_i = 0; off_i < ptElectron.size(); ++off_i) {
0187         if (sqrt((Electron->phi() - phiElectron[off_i]) * (Electron->phi() - phiElectron[off_i]) +
0188                  (Electron->eta() - etaElectron[off_i]) * (Electron->eta() - etaElectron[off_i])) < 0.5) {
0189           indexOfMatchedElectron[matchedCounter] = offlineCounter;
0190           matchedCounter++;
0191           break;
0192         }
0193       }
0194       offlineCounter++;
0195     }
0196 
0197     float pfHT = 0.0;
0198     for (reco::PFJetCollection::const_iterator i_pfjet = pfJetCollection->begin(); i_pfjet != pfJetCollection->end();
0199          ++i_pfjet) {
0200       if (i_pfjet->pt() < ptThrJet_)
0201         continue;
0202       if (fabs(i_pfjet->eta()) > etaThrJet_)
0203         continue;
0204       pfHT += i_pfjet->pt();
0205     }
0206     for (reco::CaloJetCollection::const_iterator i_calojet = caloJetCollection->begin();
0207          i_calojet != caloJetCollection->end();
0208          ++i_calojet) {
0209       if (i_calojet->pt() < ptThrJet_)
0210         continue;
0211       if (fabs(i_calojet->eta()) > etaThrJet_)
0212         continue;
0213     }
0214 
0215     if (hasFiredAuxiliaryForElectronLeg && ElectronCollection->size() > 1) {
0216       if (hasFired && indexOfMatchedElectron[1] >= 0) {  // fill trailing leg

0217         h_EleTurnOn_num->Fill(ElectronCollection->at(indexOfMatchedElectron[1]).pt());
0218       }
0219       h_EleTurnOn_den->Fill(ElectronCollection->at(1).pt());
0220     }
0221     if (hasFiredAuxiliaryForHadronicLeg) {
0222       if (hasFired) {
0223         h_pfHTTurnOn_num->Fill(pfHT);
0224       }
0225       h_pfHTTurnOn_den->Fill(pfHT);
0226     }
0227   }
0228 }
0229 
0230 void SUSY_HLT_DoubleEle_Hadronic::bookHistos(DQMStore::IBooker &ibooker_) {
0231   ibooker_.cd();
0232   ibooker_.setCurrentFolder("HLT/SUSYBSM/" + triggerPath_);
0233 
0234   // offline quantities

0235 
0236   // online quantities

0237   h_triggerElePt = ibooker_.book1D("triggerElePt", "Trigger Electron Pt; GeV", 50, 0.0, 500.0);
0238   h_triggerEleEta = ibooker_.book1D("triggerEleEta", "Trigger Electron Eta", 20, -3.0, 3.0);
0239   h_triggerElePhi = ibooker_.book1D("triggerElePhi", "Trigger Electron Phi", 20, -3.5, 3.5);
0240 
0241   h_triggerDoubleEleMass = ibooker_.book1D("triggerDoubleEleMass", "Trigger DoubleElectron Mass", 202, -2, 200);
0242 
0243   // num and den hists to be divided in harvesting step to make turn on curves

0244   h_pfHTTurnOn_num = ibooker_.book1D("pfHTTurnOn_num", "PF HT Turn On Numerator", 30, 0.0, 1500.0);
0245   h_pfHTTurnOn_den = ibooker_.book1D("pfHTTurnOn_den", "PF HT Turn On Denominator", 30, 0.0, 1500.0);
0246   h_EleTurnOn_num = ibooker_.book1D("EleTurnOn_num", "Electron Turn On Numerator", 30, 0.0, 150);
0247   h_EleTurnOn_den = ibooker_.book1D("EleTurnOn_den", "Electron Turn On Denominator", 30, 0.0, 150.0);
0248 
0249   ibooker_.cd();
0250 }
0251 
0252 // define this as a plug-in

0253 DEFINE_FWK_MODULE(SUSY_HLT_DoubleEle_Hadronic);