Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:01:14

0001 #include "FWCore/Common/interface/TriggerNames.h"
0002 #include "FWCore/Common/interface/TriggerResultsByName.h"
0003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0004 #include "DataFormats/Candidate/interface/CandidateFwd.h"
0005 #include "DataFormats/Candidate/interface/Candidate.h"
0006 #include "DataFormats/Candidate/interface/LeafCandidate.h"
0007 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
0008 
0009 #include "TString.h"
0010 #include "TRegexp.h"
0011 
0012 #include <cassert>
0013 
0014 #include "CommonTools/TriggerUtils/interface/CandidateTriggerObjectProducer.h"
0015 
0016 //
0017 // constructors and destructor
0018 //
0019 CandidateTriggerObjectProducer::CandidateTriggerObjectProducer(const edm::ParameterSet& ps)
0020     : triggerResultsToken_(consumes<edm::TriggerResults>(ps.getParameter<edm::InputTag>("triggerResults"))),
0021       triggerEventToken_(consumes<trigger::TriggerEvent>(ps.getParameter<edm::InputTag>("triggerEvent"))),
0022       processName_(ps.getParameter<edm::InputTag>("triggerEvent").process()),
0023       triggerName_(ps.getParameter<std::string>("triggerName")),
0024       hltPrescaleProvider_(ps, consumesCollector(), *this) {
0025   //   cout << "Trigger Object Producer:" << endl
0026   //        << "   TriggerResultsTag = " << triggerResultsTag_.encode() << endl
0027   //        << "   TriggerEventTag = " << triggerEventTag_.encode() << endl;
0028 
0029   produces<reco::CandidateCollection>();
0030 }
0031 
0032 CandidateTriggerObjectProducer::~CandidateTriggerObjectProducer() {}
0033 
0034 //
0035 // member functions
0036 //
0037 void CandidateTriggerObjectProducer::beginRun(const edm::Run& iRun, edm::EventSetup const& iSetup) {
0038   using namespace edm;
0039 
0040   bool changed(false);
0041   if (!hltPrescaleProvider_.init(iRun, iSetup, processName_, changed)) {
0042     edm::LogError("CandidateTriggerObjectProducer") << "Error! Can't initialize HLTConfigProvider";
0043     throw cms::Exception("HLTConfigProvider::init() returned non 0");
0044   }
0045 
0046   return;
0047 }
0048 
0049 // ------------ method called to produce the data  ------------
0050 void CandidateTriggerObjectProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0051   using namespace edm;
0052   using namespace reco;
0053   using namespace trigger;
0054 
0055   std::unique_ptr<reco::CandidateCollection> coll(new reco::CandidateCollection);
0056 
0057   // get event products
0058   edm::Handle<edm::TriggerResults> triggerResultsHandle;
0059   iEvent.getByToken(triggerResultsToken_, triggerResultsHandle);
0060   if (!triggerResultsHandle.isValid()) {
0061     edm::LogError("CandidateTriggerObjectProducer")
0062         << "CandidateTriggerObjectProducer::analyze: Error in getting TriggerResults product from Event!";
0063     return;
0064   }
0065   edm::Handle<trigger::TriggerEvent> triggerEventHandle;
0066   iEvent.getByToken(triggerEventToken_, triggerEventHandle);
0067   if (!triggerEventHandle.isValid()) {
0068     edm::LogError("CandidateTriggerObjectProducer")
0069         << "CandidateTriggerObjectProducer::analyze: Error in getting TriggerEvent product from Event!";
0070     return;
0071   }
0072 
0073   HLTConfigProvider const& hltConfig = hltPrescaleProvider_.hltConfigProvider();
0074 
0075   // sanity check
0076   //   std::cout << hltConfig.size() << std::endl;
0077   //   std::cout << triggerResultsHandle->size() << std::endl;
0078   assert(triggerResultsHandle->size() == hltConfig.size());
0079 
0080   const unsigned int n(hltConfig.size());
0081   std::vector<std::string> activeHLTPathsInThisEvent = hltConfig.triggerNames();
0082   std::map<std::string, bool> triggerInMenu;
0083   std::map<std::string, bool> triggerUnprescaled;
0084 
0085   for (std::vector<std::string>::const_iterator iHLT = activeHLTPathsInThisEvent.begin();
0086        iHLT != activeHLTPathsInThisEvent.end();
0087        ++iHLT) {
0088     //matching with regexp filter name. More than 1 matching filter is allowed
0089     if (TString(*iHLT).Contains(TRegexp(TString(triggerName_)))) {
0090       triggerInMenu[*iHLT] = true;
0091       auto const prescales = hltPrescaleProvider_.prescaleValues<double>(iEvent, iSetup, *iHLT);
0092       if (prescales.first == 1 and prescales.second == 1)
0093         triggerUnprescaled[*iHLT] = true;
0094     }
0095   }
0096 
0097   for (std::map<std::string, bool>::const_iterator iMyHLT = triggerInMenu.begin(); iMyHLT != triggerInMenu.end();
0098        ++iMyHLT) {
0099     //using only unprescaled triggers
0100     if (!(iMyHLT->second && triggerUnprescaled[iMyHLT->first]))
0101       continue;
0102     const unsigned int triggerIndex(hltConfig.triggerIndex(iMyHLT->first));
0103 
0104     assert(triggerIndex == iEvent.triggerNames(*triggerResultsHandle).triggerIndex(iMyHLT->first));
0105 
0106     // abort on invalid trigger name
0107     if (triggerIndex >= n) {
0108       edm::LogError("CandidateTriggerObjectProducer")
0109           << "CandidateTriggerObjectProducer::analyzeTrigger: path " << triggerName_ << " - not found!";
0110       return;
0111     }
0112 
0113     // modules on this trigger path
0114     //      const unsigned int m(hltConfig.size(triggerIndex));
0115     const std::vector<std::string>& moduleLabels(hltConfig.saveTagsModules(triggerIndex));
0116 
0117     // Results from TriggerResults product
0118     if (!(triggerResultsHandle->wasrun(triggerIndex)) || !(triggerResultsHandle->accept(triggerIndex)) ||
0119         (triggerResultsHandle->error(triggerIndex))) {
0120       continue;
0121     }
0122 
0123     //       const unsigned int moduleIndex(triggerResultsHandle->index(triggerIndex));
0124 
0125     //       assert (moduleIndex<m);
0126 
0127     // Results from TriggerEvent product - Looking only on last filter since trigger is accepted
0128     for (unsigned int imodule = 0; imodule < moduleLabels.size(); ++imodule) {
0129       const std::string& moduleLabel(moduleLabels[imodule]);
0130       const std::string moduleType(hltConfig.moduleType(moduleLabel));
0131       //Avoiding L1 seeds
0132       if (moduleType.find("Level1GTSeed") != std::string::npos)
0133         continue;
0134       // check whether the module is packed up in TriggerEvent product
0135       const unsigned int filterIndex(triggerEventHandle->filterIndex(InputTag(moduleLabel, "", processName_)));
0136       if (filterIndex < triggerEventHandle->sizeFilters()) {
0137         //      std::cout << " 'L3' filter in slot " << imodule
0138         //                << " - label/type " << moduleLabel << "/" << moduleType << std::endl;
0139         const Vids& VIDS(triggerEventHandle->filterIds(filterIndex));
0140         const Keys& KEYS(triggerEventHandle->filterKeys(filterIndex));
0141         const size_type nI(VIDS.size());
0142         const size_type nK(KEYS.size());
0143         assert(nI == nK);
0144         const size_type n(std::max(nI, nK));
0145         //      std::cout << "   " << n  << " accepted 'L3' objects found: " << std::endl;
0146         const TriggerObjectCollection& TOC(triggerEventHandle->getObjects());
0147         for (size_type i = 0; i != n; ++i) {
0148           const TriggerObject& TO(TOC[KEYS[i]]);
0149           coll->push_back(reco::LeafCandidate(0, TO.particle().p4(), reco::Particle::Point(0., 0., 0.), TO.id()));
0150           //          std::cout << "   " << i << " " << VIDS[i] << "/" << KEYS[i] << ": "
0151           //                      << TO.id() << " " << TO.pt() << " " << TO.eta() << " " << TO.phi() << " " << TO.mass()
0152           //                      << std::endl;
0153         }
0154       }
0155     }
0156   }
0157 
0158   iEvent.put(std::move(coll));
0159   return;
0160 }
0161 
0162 #include "FWCore/Framework/interface/MakerMacros.h"
0163 
0164 DEFINE_FWK_MODULE(CandidateTriggerObjectProducer);