Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-07-02 00:53:49

0001 /** \class HLTEgammaCombMassFilter
0002  *
0003  * $Id: HLTEgammaCombMassFilter.cc,
0004  *
0005  *  \author Chris Tully (Princeton)
0006  *
0007  */
0008 
0009 #include "HLTEgammaCombMassFilter.h"
0010 
0011 #include "DataFormats/Common/interface/Handle.h"
0012 
0013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0014 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0016 
0017 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
0018 #include "DataFormats/EgammaCandidates/interface/Electron.h"
0019 
0020 //
0021 // constructors and destructor
0022 //
0023 HLTEgammaCombMassFilter::HLTEgammaCombMassFilter(const edm::ParameterSet& iConfig) : HLTFilter(iConfig) {
0024   firstLegLastFilterTag_ = iConfig.getParameter<edm::InputTag>("firstLegLastFilter");
0025   secondLegLastFilterTag_ = iConfig.getParameter<edm::InputTag>("secondLegLastFilter");
0026   minMass_ = iConfig.getParameter<double>("minMass");
0027   l1EGTag_ = iConfig.getParameter<edm::InputTag>("l1EGCand");
0028   firstLegLastFilterToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(firstLegLastFilterTag_);
0029   secondLegLastFilterToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(secondLegLastFilterTag_);
0030 }
0031 
0032 HLTEgammaCombMassFilter::~HLTEgammaCombMassFilter() = default;
0033 
0034 void HLTEgammaCombMassFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0035   edm::ParameterSetDescription desc;
0036   makeHLTFilterDescription(desc);
0037   desc.add<edm::InputTag>("firstLegLastFilter", edm::InputTag("firstFilter"));
0038   desc.add<edm::InputTag>("secondLegLastFilter", edm::InputTag("secondFilter"));
0039   desc.add<edm::InputTag>("l1EGCand", edm::InputTag("hltEgammaCandidates"));
0040   desc.add<double>("minMass", -1.0);
0041   descriptions.add("hltEgammaCombMassFilter", desc);
0042 }
0043 
0044 // ------------ method called to produce the data  ------------
0045 bool HLTEgammaCombMassFilter::hltFilter(edm::Event& iEvent,
0046                                         const edm::EventSetup& iSetup,
0047                                         trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0048   using namespace trigger;
0049   if (saveTags()) {
0050     filterproduct.addCollectionTag(l1EGTag_);
0051   }
0052 
0053   std::vector<math::XYZTLorentzVector> firstLegP4s;
0054   std::vector<math::XYZTLorentzVector> secondLegP4s;
0055 
0056   getP4OfLegCands(iEvent, firstLegLastFilterToken_, firstLegP4s);
0057   getP4OfLegCands(iEvent, secondLegLastFilterToken_, secondLegP4s);
0058 
0059   bool accept = false;
0060   std::set<math::XYZTLorentzVector, LorentzVectorComparator> addedLegP4s;
0061 
0062   for (size_t i = 0; i < firstLegP4s.size(); i++) {
0063     for (size_t j = 0; j < secondLegP4s.size(); j++) {
0064       // Skip if it's the same object
0065       if (firstLegP4s[i] == secondLegP4s[j])
0066         continue;
0067 
0068       math::XYZTLorentzVector pairP4 = firstLegP4s[i] + secondLegP4s[j];
0069       double mass = pairP4.M();
0070       if (mass >= minMass_) {
0071         accept = true;
0072 
0073         // Add first leg object if not already added
0074         if (addedLegP4s.insert(firstLegP4s[i]).second) {
0075           addObjectToFilterProduct(iEvent, filterproduct, firstLegLastFilterToken_, i);
0076         }
0077 
0078         // Add second leg object if not already added
0079         if (addedLegP4s.insert(secondLegP4s[j]).second) {
0080           addObjectToFilterProduct(iEvent, filterproduct, secondLegLastFilterToken_, j);
0081         }
0082       }
0083     }
0084   }
0085 
0086   return accept;
0087 }
0088 
0089 void HLTEgammaCombMassFilter::addObjectToFilterProduct(
0090     const edm::Event& iEvent,
0091     trigger::TriggerFilterObjectWithRefs& filterproduct,
0092     const edm::EDGetTokenT<trigger::TriggerFilterObjectWithRefs>& token,
0093     size_t index) {
0094   edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
0095   iEvent.getByToken(token, PrevFilterOutput);
0096 
0097   // Get all types of objects
0098   std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > phoCandsPrev;
0099   PrevFilterOutput->getObjects(trigger::TriggerPhoton, phoCandsPrev);
0100   std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > clusCandsPrev;
0101   PrevFilterOutput->getObjects(trigger::TriggerCluster, clusCandsPrev);
0102   std::vector<edm::Ref<reco::ElectronCollection> > eleCandsPrev;
0103   PrevFilterOutput->getObjects(trigger::TriggerElectron, eleCandsPrev);
0104 
0105   // Check which type of object corresponds to the given index
0106   if (index < phoCandsPrev.size()) {
0107     filterproduct.addObject(trigger::TriggerPhoton, phoCandsPrev[index]);
0108   } else if (index < phoCandsPrev.size() + clusCandsPrev.size()) {
0109     filterproduct.addObject(trigger::TriggerCluster, clusCandsPrev[index - phoCandsPrev.size()]);
0110   } else if (index < phoCandsPrev.size() + clusCandsPrev.size() + eleCandsPrev.size()) {
0111     filterproduct.addObject(trigger::TriggerElectron, eleCandsPrev[index - phoCandsPrev.size() - clusCandsPrev.size()]);
0112   } else {
0113     edm::LogWarning("HLTEgammaCombMassFilter") << "Could not find object at index " << index;
0114   }
0115 }
0116 
0117 void HLTEgammaCombMassFilter::getP4OfLegCands(const edm::Event& iEvent,
0118                                               const edm::EDGetTokenT<trigger::TriggerFilterObjectWithRefs>& filterToken,
0119                                               std::vector<math::XYZTLorentzVector>& p4s) {
0120   edm::Handle<trigger::TriggerFilterObjectWithRefs> filterOutput;
0121   iEvent.getByToken(filterToken, filterOutput);
0122 
0123   //its easier on the if statement flow if I try everything at once, shouldnt add to timing
0124   std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > phoCands;
0125   filterOutput->getObjects(trigger::TriggerPhoton, phoCands);
0126   std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > clusCands;
0127   filterOutput->getObjects(trigger::TriggerCluster, clusCands);
0128   std::vector<edm::Ref<reco::ElectronCollection> > eleCands;
0129   filterOutput->getObjects(trigger::TriggerElectron, eleCands);
0130 
0131   if (!phoCands.empty()) {  //its photons
0132     for (auto const& phoCand : phoCands) {
0133       p4s.push_back(phoCand->p4());
0134     }
0135   } else if (!clusCands.empty()) {
0136     //try trigger cluster (should never be this, at the time of writing (17/1/11) this would indicate an error)
0137     for (auto const& clusCand : clusCands) {
0138       p4s.push_back(clusCand->p4());
0139     }
0140   } else if (!eleCands.empty()) {
0141     for (auto const& eleCand : eleCands) {
0142       p4s.push_back(eleCand->p4());
0143     }
0144   }
0145 }
0146 
0147 // declare this class as a framework plugin
0148 #include "FWCore/Framework/interface/MakerMacros.h"
0149 DEFINE_FWK_MODULE(HLTEgammaCombMassFilter);