Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2022-02-16 06:15:45

0001 /** \class HLTEventAnalyzerRAW
0002  *
0003  * See header file for documentation
0004  *
0005  *
0006  *  \author Martin Grunewald
0007  *
0008  */
0009 
0010 #include "FWCore/Common/interface/TriggerNames.h"
0011 #include "FWCore/Common/interface/TriggerResultsByName.h"
0012 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0013 #include "HLTrigger/HLTcore/interface/HLTEventAnalyzerRAW.h"
0014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0015 
0016 // need access to class objects being referenced to get their content!
0017 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
0018 #include "DataFormats/EgammaCandidates/interface/Electron.h"
0019 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
0020 #include "DataFormats/JetReco/interface/CaloJet.h"
0021 #include "DataFormats/Candidate/interface/CompositeCandidate.h"
0022 #include "DataFormats/METReco/interface/MET.h"
0023 #include "DataFormats/METReco/interface/CaloMET.h"
0024 #include "DataFormats/METReco/interface/PFMET.h"
0025 #include "DataFormats/HcalIsolatedTrack/interface/IsolatedPixelTrackCandidate.h"
0026 
0027 #include "DataFormats/L1Trigger/interface/L1HFRings.h"         // deprecate
0028 #include "DataFormats/L1Trigger/interface/L1EmParticle.h"      // deprecate
0029 #include "DataFormats/L1Trigger/interface/L1JetParticle.h"     // deprecate
0030 #include "DataFormats/L1Trigger/interface/L1MuonParticle.h"    // deprecate
0031 #include "DataFormats/L1Trigger/interface/L1EtMissParticle.h"  // deprecate
0032 
0033 #include "DataFormats/L1Trigger/interface/Muon.h"
0034 #include "DataFormats/L1Trigger/interface/MuonShower.h"
0035 #include "DataFormats/L1Trigger/interface/EGamma.h"
0036 #include "DataFormats/L1Trigger/interface/Jet.h"
0037 #include "DataFormats/L1Trigger/interface/Tau.h"
0038 #include "DataFormats/L1Trigger/interface/EtSum.h"
0039 
0040 #include "DataFormats/L1TCorrelator/interface/TkMuon.h"
0041 #include "DataFormats/L1TCorrelator/interface/TkElectron.h"
0042 #include "DataFormats/L1TCorrelator/interface/TkEm.h"
0043 #include "DataFormats/L1TParticleFlow/interface/PFJet.h"
0044 #include "DataFormats/L1TParticleFlow/interface/PFTau.h"
0045 #include "DataFormats/L1TParticleFlow/interface/HPSPFTau.h"
0046 #include "DataFormats/L1TParticleFlow/interface/PFTrack.h"
0047 #include "DataFormats/JetReco/interface/PFJet.h"
0048 #include "DataFormats/TauReco/interface/PFTau.h"
0049 
0050 #include <cassert>
0051 
0052 //
0053 // constructors and destructor
0054 //
0055 HLTEventAnalyzerRAW::HLTEventAnalyzerRAW(const edm::ParameterSet& ps)
0056     : processName_(ps.getParameter<std::string>("processName")),
0057       triggerName_(ps.getParameter<std::string>("triggerName")),
0058       triggerResultsTag_(ps.getParameter<edm::InputTag>("triggerResults")),
0059       triggerResultsToken_(consumes<edm::TriggerResults>(triggerResultsTag_)),
0060       triggerEventWithRefsTag_(ps.getParameter<edm::InputTag>("triggerEventWithRefs")),
0061       triggerEventWithRefsToken_(consumes<trigger::TriggerEventWithRefs>(triggerEventWithRefsTag_)) {
0062   using namespace std;
0063   using namespace edm;
0064 
0065   LogVerbatim("HLTEventAnalyzerRAW") << "HLTEventAnalyzerRAW configuration: " << endl
0066                                      << "   ProcessName = " << processName_ << endl
0067                                      << "   TriggerName = " << triggerName_ << endl
0068                                      << "   TriggerResultsTag = " << triggerResultsTag_.encode() << endl
0069                                      << "   TriggerEventWithRefsTag = " << triggerEventWithRefsTag_.encode() << endl;
0070 }
0071 
0072 HLTEventAnalyzerRAW::~HLTEventAnalyzerRAW() = default;
0073 
0074 //
0075 // member functions
0076 //
0077 void HLTEventAnalyzerRAW::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0078   edm::ParameterSetDescription desc;
0079   desc.add<std::string>("processName", "HLT");
0080   desc.add<std::string>("triggerName", "@");
0081   desc.add<edm::InputTag>("triggerResults", edm::InputTag("TriggerResults", "", "HLT"));
0082   desc.add<edm::InputTag>("triggerEventWithRefs", edm::InputTag("hltTriggerSummaryRAW", "", "HLT"));
0083   descriptions.add("hltEventAnalyzerRAW", desc);
0084 }
0085 
0086 void HLTEventAnalyzerRAW::endRun(edm::Run const& iRun, edm::EventSetup const& iSetup) {}
0087 
0088 void HLTEventAnalyzerRAW::beginRun(edm::Run const& iRun, edm::EventSetup const& iSetup) {
0089   using namespace std;
0090   using namespace edm;
0091 
0092   bool changed(true);
0093   if (hltConfig_.init(iRun, iSetup, processName_, changed)) {
0094     if (changed) {
0095       // check if trigger name in (new) config
0096       if (triggerName_ != "@") {  // "@" means: analyze all triggers in config
0097         const unsigned int n(hltConfig_.size());
0098         const unsigned int triggerIndex(hltConfig_.triggerIndex(triggerName_));
0099         if (triggerIndex >= n) {
0100           LogVerbatim("HLTEventAnalyzerRAW")
0101               << "HLTEventAnalyzerRAW::analyze:"
0102               << " TriggerName " << triggerName_ << " not available in (new) config!" << endl;
0103           LogVerbatim("HLTEventAnalyzerRAW") << "Available TriggerNames are: " << endl;
0104           hltConfig_.dump("Triggers");
0105         }
0106       }
0107     }
0108   } else {
0109     LogVerbatim("HLTEventAnalyzerRAW") << "HLTEventAnalyzerRAW::analyze:"
0110                                        << " config extraction failure with process name " << processName_ << endl;
0111   }
0112 }
0113 
0114 // ------------ method called to produce the data  ------------
0115 void HLTEventAnalyzerRAW::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0116   using namespace std;
0117   using namespace edm;
0118 
0119   LogVerbatim("HLTEventAnalyzerRAW") << endl;
0120 
0121   // get event products
0122   iEvent.getByToken(triggerResultsToken_, triggerResultsHandle_);
0123   if (!triggerResultsHandle_.isValid()) {
0124     LogVerbatim("HLTEventAnalyzerRAW")
0125         << "HLTEventAnalyzerRAW::analyze: Error in getting TriggerResults product from Event!" << endl;
0126     return;
0127   }
0128   iEvent.getByToken(triggerEventWithRefsToken_, triggerEventWithRefsHandle_);
0129   if (!triggerEventWithRefsHandle_.isValid()) {
0130     LogVerbatim("HLTEventAnalyzerRAW")
0131         << "HLTEventAnalyzerRAW::analyze: Error in getting TriggerEventWithRefs product from Event!" << endl;
0132     return;
0133   }
0134   // sanity check
0135   assert(triggerResultsHandle_->size() == hltConfig_.size());
0136 
0137   // analyze this event for the triggers requested
0138   if (triggerName_ == "@") {
0139     const unsigned int n(hltConfig_.size());
0140     for (unsigned int i = 0; i != n; ++i) {
0141       analyzeTrigger(iEvent, iSetup, hltConfig_.triggerName(i));
0142     }
0143   } else {
0144     analyzeTrigger(iEvent, iSetup, triggerName_);
0145   }
0146 
0147   return;
0148 }
0149 
0150 void HLTEventAnalyzerRAW::analyzeTrigger(const edm::Event& iEvent,
0151                                          const edm::EventSetup& iSetup,
0152                                          const std::string& triggerName) {
0153   using namespace std;
0154   using namespace edm;
0155   using namespace reco;
0156   using namespace trigger;
0157 
0158   LogVerbatim("HLTEventAnalyzerRAW") << endl;
0159 
0160   const unsigned int n(hltConfig_.size());
0161   const unsigned int triggerIndex(hltConfig_.triggerIndex(triggerName));
0162   assert(triggerIndex == iEvent.triggerNames(*triggerResultsHandle_).triggerIndex(triggerName));
0163 
0164   // abort on invalid trigger name
0165   if (triggerIndex >= n) {
0166     LogVerbatim("HLTEventAnalyzerRAW") << "HLTEventAnalyzerRAW::analyzeTrigger: path " << triggerName << " - not found!"
0167                                        << endl;
0168     return;
0169   }
0170 
0171   LogVerbatim("HLTEventAnalyzerRAW") << "HLTEventAnalyzerRAW::analyzeTrigger: path " << triggerName << " ["
0172                                      << triggerIndex << "]" << endl;
0173   // modules on this trigger path
0174   const unsigned int m(hltConfig_.size(triggerIndex));
0175   const vector<string>& moduleLabels(hltConfig_.moduleLabels(triggerIndex));
0176 
0177   // Results from TriggerResults product
0178   LogVerbatim("HLTEventAnalyzerRAW") << " Trigger path status:"
0179                                      << " WasRun=" << triggerResultsHandle_->wasrun(triggerIndex)
0180                                      << " Accept=" << triggerResultsHandle_->accept(triggerIndex)
0181                                      << " Error =" << triggerResultsHandle_->error(triggerIndex) << endl;
0182   const unsigned int moduleIndex(triggerResultsHandle_->index(triggerIndex));
0183   LogVerbatim("HLTEventAnalyzerRAW") << " Last active module - label/type: " << moduleLabels[moduleIndex] << "/"
0184                                      << hltConfig_.moduleType(moduleLabels[moduleIndex]) << " [" << moduleIndex
0185                                      << " out of 0-" << (m - 1) << " on this path]" << endl;
0186   assert(moduleIndex < m);
0187 
0188   // Results from TriggerEventWithRefs product
0189   photonIds_.clear();
0190   photonRefs_.clear();
0191   electronIds_.clear();
0192   electronRefs_.clear();
0193   muonIds_.clear();
0194   muonRefs_.clear();
0195   jetIds_.clear();
0196   jetRefs_.clear();
0197   compositeIds_.clear();
0198   compositeRefs_.clear();
0199   basemetIds_.clear();
0200   basemetRefs_.clear();
0201   calometIds_.clear();
0202   calometRefs_.clear();
0203   pixtrackIds_.clear();
0204   pixtrackRefs_.clear();
0205 
0206   l1emIds_.clear();
0207   l1emRefs_.clear();
0208   l1muonIds_.clear();
0209   l1muonRefs_.clear();
0210   l1jetIds_.clear();
0211   l1jetRefs_.clear();
0212   l1etmissIds_.clear();
0213   l1etmissRefs_.clear();
0214   l1hfringsIds_.clear();
0215   l1hfringsRefs_.clear();
0216 
0217   l1tmuonIds_.clear();
0218   l1tmuonRefs_.clear();
0219   l1tmuonShowerIds_.clear();
0220   l1tmuonShowerRefs_.clear();
0221   l1tegammaIds_.clear();
0222   l1tegammaRefs_.clear();
0223   l1tjetIds_.clear();
0224   l1tjetRefs_.clear();
0225   l1ttauIds_.clear();
0226   l1ttauRefs_.clear();
0227   l1tetsumIds_.clear();
0228   l1tetsumRefs_.clear();
0229 
0230   l1ttkmuIds_.clear();
0231   l1ttkmuRefs_.clear();
0232   l1ttkeleIds_.clear();
0233   l1ttkeleRefs_.clear();
0234   l1ttkemIds_.clear();
0235   l1ttkemRefs_.clear();
0236   l1tpfjetIds_.clear();
0237   l1tpfjetRefs_.clear();
0238   l1tpftauIds_.clear();
0239   l1tpftauRefs_.clear();
0240   l1thpspftauIds_.clear();
0241   l1thpspftauRefs_.clear();
0242   l1tpftrackIds_.clear();
0243   l1tpftrackRefs_.clear();
0244 
0245   pfjetIds_.clear();
0246   pfjetRefs_.clear();
0247   pftauIds_.clear();
0248   pftauRefs_.clear();
0249   pfmetIds_.clear();
0250   pfmetRefs_.clear();
0251 
0252   // Attention: must look only for modules actually run in this path
0253   // for this event!
0254   for (unsigned int j = 0; j <= moduleIndex; ++j) {
0255     const string& moduleLabel(moduleLabels[j]);
0256     const string moduleType(hltConfig_.moduleType(moduleLabel));
0257     // check whether the module is packed up in TriggerEventWithRef product
0258     const unsigned int filterIndex(triggerEventWithRefsHandle_->filterIndex(InputTag(moduleLabel, "", processName_)));
0259     if (filterIndex < triggerEventWithRefsHandle_->size()) {
0260       LogVerbatim("HLTEventAnalyzerRAW") << " Filter in slot " << j << " - label/type " << moduleLabel << "/"
0261                                          << moduleType << endl;
0262       LogVerbatim("HLTEventAnalyzerRAW") << " Filter packed up at: " << filterIndex << endl;
0263       LogVerbatim("HLTEventAnalyzerRAW") << "  Accepted objects:" << endl;
0264 
0265       triggerEventWithRefsHandle_->getObjects(filterIndex, photonIds_, photonRefs_);
0266       const unsigned int nPhotons(photonIds_.size());
0267       if (nPhotons > 0) {
0268         LogVerbatim("HLTEventAnalyzerRAW") << "   Photons: " << nPhotons << "  - the objects: # id pt" << endl;
0269         for (unsigned int i = 0; i != nPhotons; ++i) {
0270           LogVerbatim("HLTEventAnalyzerRAW")
0271               << "   " << i << " " << photonIds_[i] << " " << photonRefs_[i]->pt() << endl;
0272         }
0273       }
0274 
0275       triggerEventWithRefsHandle_->getObjects(filterIndex, electronIds_, electronRefs_);
0276       const unsigned int nElectrons(electronIds_.size());
0277       if (nElectrons > 0) {
0278         LogVerbatim("HLTEventAnalyzerRAW") << "   Electrons: " << nElectrons << "  - the objects: # id pt" << endl;
0279         for (unsigned int i = 0; i != nElectrons; ++i) {
0280           LogVerbatim("HLTEventAnalyzerRAW")
0281               << "   " << i << " " << electronIds_[i] << " " << electronRefs_[i]->pt() << endl;
0282         }
0283       }
0284 
0285       triggerEventWithRefsHandle_->getObjects(filterIndex, muonIds_, muonRefs_);
0286       const unsigned int nMuons(muonIds_.size());
0287       if (nMuons > 0) {
0288         LogVerbatim("HLTEventAnalyzerRAW") << "   Muons: " << nMuons << "  - the objects: # id pt" << endl;
0289         for (unsigned int i = 0; i != nMuons; ++i) {
0290           LogVerbatim("HLTEventAnalyzerRAW") << "   " << i << " " << muonIds_[i] << " " << muonRefs_[i]->pt() << endl;
0291         }
0292       }
0293 
0294       triggerEventWithRefsHandle_->getObjects(filterIndex, jetIds_, jetRefs_);
0295       const unsigned int nJets(jetIds_.size());
0296       if (nJets > 0) {
0297         LogVerbatim("HLTEventAnalyzerRAW") << "   Jets: " << nJets << "  - the objects: # id pt" << endl;
0298         for (unsigned int i = 0; i != nJets; ++i) {
0299           LogVerbatim("HLTEventAnalyzerRAW") << "   " << i << " " << jetIds_[i] << " " << jetRefs_[i]->pt() << endl;
0300         }
0301       }
0302 
0303       triggerEventWithRefsHandle_->getObjects(filterIndex, compositeIds_, compositeRefs_);
0304       const unsigned int nComposites(compositeIds_.size());
0305       if (nComposites > 0) {
0306         LogVerbatim("HLTEventAnalyzerRAW") << "   Composites: " << nComposites << "  - the objects: # id pt" << endl;
0307         for (unsigned int i = 0; i != nComposites; ++i) {
0308           LogVerbatim("HLTEventAnalyzerRAW")
0309               << "   " << i << " " << compositeIds_[i] << " " << compositeRefs_[i]->pt() << endl;
0310         }
0311       }
0312 
0313       triggerEventWithRefsHandle_->getObjects(filterIndex, basemetIds_, basemetRefs_);
0314       const unsigned int nBaseMETs(basemetIds_.size());
0315       if (nBaseMETs > 0) {
0316         LogVerbatim("HLTEventAnalyzerRAW") << "   BaseMETs: " << nBaseMETs << "  - the objects: # id pt" << endl;
0317         for (unsigned int i = 0; i != nBaseMETs; ++i) {
0318           LogVerbatim("HLTEventAnalyzerRAW")
0319               << "   " << i << " " << basemetIds_[i] << " " << basemetRefs_[i]->pt() << endl;
0320         }
0321       }
0322 
0323       triggerEventWithRefsHandle_->getObjects(filterIndex, calometIds_, calometRefs_);
0324       const unsigned int nCaloMETs(calometIds_.size());
0325       if (nCaloMETs > 0) {
0326         LogVerbatim("HLTEventAnalyzerRAW") << "   CaloMETs: " << nCaloMETs << "  - the objects: # id pt" << endl;
0327         for (unsigned int i = 0; i != nCaloMETs; ++i) {
0328           LogVerbatim("HLTEventAnalyzerRAW")
0329               << "   " << i << " " << calometIds_[i] << " " << calometRefs_[i]->pt() << endl;
0330         }
0331       }
0332 
0333       triggerEventWithRefsHandle_->getObjects(filterIndex, pixtrackIds_, pixtrackRefs_);
0334       const unsigned int nPixTracks(pixtrackIds_.size());
0335       if (nPixTracks > 0) {
0336         LogVerbatim("HLTEventAnalyzerRAW") << "   PixTracks: " << nPixTracks << "  - the objects: # id pt" << endl;
0337         for (unsigned int i = 0; i != nPixTracks; ++i) {
0338           LogVerbatim("HLTEventAnalyzerRAW")
0339               << "   " << i << " " << pixtrackIds_[i] << " " << pixtrackRefs_[i]->pt() << endl;
0340         }
0341       }
0342 
0343       triggerEventWithRefsHandle_->getObjects(filterIndex, l1emIds_, l1emRefs_);
0344       const unsigned int nL1EM(l1emIds_.size());
0345       if (nL1EM > 0) {
0346         LogVerbatim("HLTEventAnalyzerRAW") << "   L1EM: " << nL1EM << "  - the objects: # id pt" << endl;
0347         for (unsigned int i = 0; i != nL1EM; ++i) {
0348           LogVerbatim("HLTEventAnalyzerRAW") << "   " << i << " " << l1emIds_[i] << " " << l1emRefs_[i]->pt() << endl;
0349         }
0350       }
0351 
0352       triggerEventWithRefsHandle_->getObjects(filterIndex, l1muonIds_, l1muonRefs_);
0353       const unsigned int nL1Muon(l1muonIds_.size());
0354       if (nL1Muon > 0) {
0355         LogVerbatim("HLTEventAnalyzerRAW") << "   L1Muon: " << nL1Muon << "  - the objects: # id pt" << endl;
0356         for (unsigned int i = 0; i != nL1Muon; ++i) {
0357           LogVerbatim("HLTEventAnalyzerRAW")
0358               << "   " << i << " " << l1muonIds_[i] << " " << l1muonRefs_[i]->pt() << endl;
0359         }
0360       }
0361 
0362       triggerEventWithRefsHandle_->getObjects(filterIndex, l1jetIds_, l1jetRefs_);
0363       const unsigned int nL1Jet(l1jetIds_.size());
0364       if (nL1Jet > 0) {
0365         LogVerbatim("HLTEventAnalyzerRAW") << "   L1Jet: " << nL1Jet << "  - the objects: # id pt" << endl;
0366         for (unsigned int i = 0; i != nL1Jet; ++i) {
0367           LogVerbatim("HLTEventAnalyzerRAW") << "   " << i << " " << l1jetIds_[i] << " " << l1jetRefs_[i]->pt() << endl;
0368         }
0369       }
0370 
0371       triggerEventWithRefsHandle_->getObjects(filterIndex, l1etmissIds_, l1etmissRefs_);
0372       const unsigned int nL1EtMiss(l1etmissIds_.size());
0373       if (nL1EtMiss > 0) {
0374         LogVerbatim("HLTEventAnalyzerRAW") << "   L1EtMiss: " << nL1EtMiss << "  - the objects: # id pt" << endl;
0375         for (unsigned int i = 0; i != nL1EtMiss; ++i) {
0376           LogVerbatim("HLTEventAnalyzerRAW")
0377               << "   " << i << " " << l1etmissIds_[i] << " " << l1etmissRefs_[i]->pt() << endl;
0378         }
0379       }
0380 
0381       triggerEventWithRefsHandle_->getObjects(filterIndex, l1hfringsIds_, l1hfringsRefs_);
0382       const unsigned int nL1HfRings(l1hfringsIds_.size());
0383       if (nL1HfRings > 0) {
0384         LogVerbatim("HLTEventAnalyzerRAW") << "   L1HfRings: " << nL1HfRings << "  - the objects: # id 4 4" << endl;
0385         for (unsigned int i = 0; i != nL1HfRings; ++i) {
0386           LogVerbatim("HLTEventAnalyzerRAW") << "   " << i << " " << l1hfringsIds_[i] << " "
0387                                              << l1hfringsRefs_[i]->hfEtSum(l1extra::L1HFRings::kRing1PosEta) << " "
0388                                              << l1hfringsRefs_[i]->hfEtSum(l1extra::L1HFRings::kRing1NegEta) << " "
0389                                              << l1hfringsRefs_[i]->hfEtSum(l1extra::L1HFRings::kRing2PosEta) << " "
0390                                              << l1hfringsRefs_[i]->hfEtSum(l1extra::L1HFRings::kRing2NegEta) << " "
0391                                              << l1hfringsRefs_[i]->hfBitCount(l1extra::L1HFRings::kRing1PosEta) << " "
0392                                              << l1hfringsRefs_[i]->hfBitCount(l1extra::L1HFRings::kRing1NegEta) << " "
0393                                              << l1hfringsRefs_[i]->hfBitCount(l1extra::L1HFRings::kRing2PosEta) << " "
0394                                              << l1hfringsRefs_[i]->hfBitCount(l1extra::L1HFRings::kRing2NegEta) << endl;
0395         }
0396       }
0397 
0398       triggerEventWithRefsHandle_->getObjects(filterIndex, l1tmuonIds_, l1tmuonRefs_);
0399       const unsigned int nL1TMuon(l1tmuonIds_.size());
0400       if (nL1TMuon > 0) {
0401         LogVerbatim("HLTEventAnalyzerRAW") << "   L1TMuon: " << nL1TMuon << "  - the objects: # id pt" << endl;
0402         for (unsigned int i = 0; i != nL1TMuon; ++i) {
0403           LogVerbatim("HLTEventAnalyzerRAW")
0404               << "   " << i << " " << l1tmuonIds_[i] << " " << l1tmuonRefs_[i]->pt() << endl;
0405         }
0406       }
0407 
0408       triggerEventWithRefsHandle_->getObjects(filterIndex, l1tmuonShowerIds_, l1tmuonShowerRefs_);
0409       const unsigned int nL1TMuonShower(l1tmuonShowerIds_.size());
0410       if (nL1TMuonShower > 0) {
0411         LogVerbatim("HLTEventAnalyzerRAW")
0412             << "   L1TMuonShower: " << nL1TMuonShower << "  - the objects: # id pt" << endl;
0413         for (unsigned int i = 0; i != nL1TMuonShower; ++i) {
0414           LogVerbatim("HLTEventAnalyzerRAW")
0415               << "   " << i << " " << l1tmuonShowerIds_[i] << " " << l1tmuonShowerRefs_[i]->pt() << endl;
0416         }
0417       }
0418 
0419       triggerEventWithRefsHandle_->getObjects(filterIndex, l1tegammaIds_, l1tegammaRefs_);
0420       const unsigned int nL1TEGamma(l1tegammaIds_.size());
0421       if (nL1TEGamma > 0) {
0422         LogVerbatim("HLTEventAnalyzerRAW") << "   L1TEGamma: " << nL1TEGamma << "  - the objects: # id pt" << endl;
0423         for (unsigned int i = 0; i != nL1TEGamma; ++i) {
0424           LogVerbatim("HLTEventAnalyzerRAW")
0425               << "   " << i << " " << l1tegammaIds_[i] << " " << l1tegammaRefs_[i]->pt() << endl;
0426         }
0427       }
0428 
0429       triggerEventWithRefsHandle_->getObjects(filterIndex, l1tjetIds_, l1tjetRefs_);
0430       const unsigned int nL1TJet(l1tjetIds_.size());
0431       if (nL1TJet > 0) {
0432         LogVerbatim("HLTEventAnalyzerRAW") << "   L1TJet: " << nL1TJet << "  - the objects: # id pt" << endl;
0433         for (unsigned int i = 0; i != nL1TJet; ++i) {
0434           LogVerbatim("HLTEventAnalyzerRAW")
0435               << "   " << i << " " << l1tjetIds_[i] << " " << l1tjetRefs_[i]->pt() << endl;
0436         }
0437       }
0438 
0439       triggerEventWithRefsHandle_->getObjects(filterIndex, l1ttauIds_, l1ttauRefs_);
0440       const unsigned int nL1TTau(l1ttauIds_.size());
0441       if (nL1TTau > 0) {
0442         LogVerbatim("HLTEventAnalyzerRAW") << "   L1TTau: " << nL1TTau << "  - the objects: # id pt" << endl;
0443         for (unsigned int i = 0; i != nL1TTau; ++i) {
0444           LogVerbatim("HLTEventAnalyzerRAW")
0445               << "   " << i << " " << l1ttauIds_[i] << " " << l1ttauRefs_[i]->pt() << endl;
0446         }
0447       }
0448 
0449       triggerEventWithRefsHandle_->getObjects(filterIndex, l1tetsumIds_, l1tetsumRefs_);
0450       const unsigned int nL1TEtSum(l1tetsumIds_.size());
0451       if (nL1TEtSum > 0) {
0452         LogVerbatim("HLTEventAnalyzerRAW") << "   L1TEtSum: " << nL1TEtSum << "  - the objects: # id pt" << endl;
0453         for (unsigned int i = 0; i != nL1TEtSum; ++i) {
0454           LogVerbatim("HLTEventAnalyzerRAW")
0455               << "   " << i << " " << l1tetsumIds_[i] << " " << l1tetsumRefs_[i]->pt() << endl;
0456         }
0457       }
0458 
0459       /* Phase-2 */
0460 
0461       triggerEventWithRefsHandle_->getObjects(filterIndex, l1ttkmuIds_, l1ttkmuRefs_);
0462       const unsigned int nL1TTkMuons(l1ttkmuIds_.size());
0463       if (nL1TTkMuons > 0) {
0464         LogVerbatim("HLTEventAnalyzerRAW") << "   L1TTkMuons: " << nL1TTkMuons << "  - the objects: # id pt" << endl;
0465         for (unsigned int i = 0; i != nL1TTkMuons; ++i) {
0466           LogVerbatim("HLTEventAnalyzerRAW")
0467               << "   " << i << " " << l1ttkmuIds_[i] << " " << l1ttkmuRefs_[i]->pt() << endl;
0468         }
0469       }
0470 
0471       triggerEventWithRefsHandle_->getObjects(filterIndex, l1ttkeleIds_, l1ttkeleRefs_);
0472       const unsigned int nL1TTkElectrons(l1ttkeleIds_.size());
0473       if (nL1TTkElectrons > 0) {
0474         LogVerbatim("HLTEventAnalyzerRAW")
0475             << "   L1TTkElectrons: " << nL1TTkElectrons << "  - the objects: # id pt" << endl;
0476         for (unsigned int i = 0; i != nL1TTkElectrons; ++i) {
0477           LogVerbatim("HLTEventAnalyzerRAW")
0478               << "   " << i << " " << l1ttkeleIds_[i] << " " << l1ttkeleRefs_[i]->pt() << endl;
0479         }
0480       }
0481 
0482       triggerEventWithRefsHandle_->getObjects(filterIndex, l1ttkemIds_, l1ttkemRefs_);
0483       const unsigned int nL1TTkEMs(l1ttkemIds_.size());
0484       if (nL1TTkEMs > 0) {
0485         LogVerbatim("HLTEventAnalyzerRAW") << "   L1TTkEMs: " << nL1TTkEMs << "  - the objects: # id pt" << endl;
0486         for (unsigned int i = 0; i != nL1TTkEMs; ++i) {
0487           LogVerbatim("HLTEventAnalyzerRAW")
0488               << "   " << i << " " << l1ttkemIds_[i] << " " << l1ttkemRefs_[i]->pt() << endl;
0489         }
0490       }
0491 
0492       triggerEventWithRefsHandle_->getObjects(filterIndex, l1tpfjetIds_, l1tpfjetRefs_);
0493       const unsigned int nL1TPFJets(l1tpfjetIds_.size());
0494       if (nL1TPFJets > 0) {
0495         LogVerbatim("HLTEventAnalyzerRAW") << "   L1TPFJets: " << nL1TPFJets << "  - the objects: # id pt" << endl;
0496         for (unsigned int i = 0; i != nL1TPFJets; ++i) {
0497           LogVerbatim("HLTEventAnalyzerRAW")
0498               << "   " << i << " " << l1tpfjetIds_[i] << " " << l1tpfjetRefs_[i]->pt() << endl;
0499         }
0500       }
0501 
0502       triggerEventWithRefsHandle_->getObjects(filterIndex, l1tpftauIds_, l1tpftauRefs_);
0503       const unsigned int nL1TPFTaus(l1tpftauIds_.size());
0504       if (nL1TPFTaus > 0) {
0505         LogVerbatim("HLTEventAnalyzerRAW") << "   L1TPFTaus: " << nL1TPFTaus << "  - the objects: # id pt" << endl;
0506         for (unsigned int i = 0; i != nL1TPFTaus; ++i) {
0507           LogVerbatim("HLTEventAnalyzerRAW")
0508               << "   " << i << " " << l1tpftauIds_[i] << " " << l1tpftauRefs_[i]->pt() << endl;
0509         }
0510       }
0511 
0512       triggerEventWithRefsHandle_->getObjects(filterIndex, l1thpspftauIds_, l1thpspftauRefs_);
0513       const unsigned int nL1THPSPFTaus(l1thpspftauIds_.size());
0514       if (nL1THPSPFTaus > 0) {
0515         LogVerbatim("HLTEventAnalyzerRAW")
0516             << "   L1THPSPFTaus: " << nL1THPSPFTaus << "  - the objects: # id pt" << endl;
0517         for (unsigned int i = 0; i != nL1THPSPFTaus; ++i) {
0518           LogVerbatim("HLTEventAnalyzerRAW")
0519               << "   " << i << " " << l1thpspftauIds_[i] << " " << l1thpspftauRefs_[i]->pt() << endl;
0520         }
0521       }
0522 
0523       triggerEventWithRefsHandle_->getObjects(filterIndex, l1tpftrackIds_, l1tpftrackRefs_);
0524       const unsigned int nL1TPFTracks(l1tpftrackIds_.size());
0525       if (nL1TPFTracks > 0) {
0526         LogVerbatim("HLTEventAnalyzerRAW") << "   L1TPFTracks: " << nL1TPFTracks << "  - the objects: # id pt" << endl;
0527         for (unsigned int i = 0; i != nL1TPFTracks; ++i) {
0528           LogVerbatim("HLTEventAnalyzerRAW")
0529               << "   " << i << " " << l1tpftrackIds_[i] << " " << l1tpftrackRefs_[i]->pt() << endl;
0530         }
0531       }
0532 
0533       triggerEventWithRefsHandle_->getObjects(filterIndex, pfjetIds_, pfjetRefs_);
0534       const unsigned int nPFJets(pfjetIds_.size());
0535       if (nPFJets > 0) {
0536         LogVerbatim("HLTEventAnalyzerRAW") << "   PFJets: " << nPFJets << "  - the objects: # id pt" << endl;
0537         for (unsigned int i = 0; i != nPFJets; ++i) {
0538           LogVerbatim("HLTEventAnalyzerRAW") << "   " << i << " " << pfjetIds_[i] << " " << pfjetRefs_[i]->pt() << endl;
0539         }
0540       }
0541 
0542       triggerEventWithRefsHandle_->getObjects(filterIndex, pftauIds_, pftauRefs_);
0543       const unsigned int nPFTaus(pftauIds_.size());
0544       if (nPFTaus > 0) {
0545         LogVerbatim("HLTEventAnalyzerRAW") << "   PFTaus: " << nPFTaus << "  - the objects: # id pt" << endl;
0546         for (unsigned int i = 0; i != nPFTaus; ++i) {
0547           LogVerbatim("HLTEventAnalyzerRAW") << "   " << i << " " << pftauIds_[i] << " " << pftauRefs_[i]->pt() << endl;
0548         }
0549       }
0550 
0551       triggerEventWithRefsHandle_->getObjects(filterIndex, pfmetIds_, pfmetRefs_);
0552       const unsigned int nPfMETs(pfmetIds_.size());
0553       if (nPfMETs > 0) {
0554         LogVerbatim("HLTEventAnalyzerRAW") << "   PfMETs: " << nPfMETs << "  - the objects: # id pt" << endl;
0555         for (unsigned int i = 0; i != nPfMETs; ++i) {
0556           LogVerbatim("HLTEventAnalyzerRAW") << "   " << i << " " << pfmetIds_[i] << " " << pfmetRefs_[i]->pt() << endl;
0557         }
0558       }
0559     }
0560   }
0561 
0562   return;
0563 }