Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:09:16

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   firstLegLastFilterToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(firstLegLastFilterTag_);
0028   secondLegLastFilterToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(secondLegLastFilterTag_);
0029 }
0030 
0031 HLTEgammaCombMassFilter::~HLTEgammaCombMassFilter() = default;
0032 
0033 void HLTEgammaCombMassFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0034   edm::ParameterSetDescription desc;
0035   makeHLTFilterDescription(desc);
0036   desc.add<edm::InputTag>("firstLegLastFilter", edm::InputTag("firstFilter"));
0037   desc.add<edm::InputTag>("secondLegLastFilter", edm::InputTag("secondFilter"));
0038   desc.add<double>("minMass", -1.0);
0039   descriptions.add("hltEgammaCombMassFilter", desc);
0040 }
0041 
0042 // ------------ method called to produce the data  ------------
0043 bool HLTEgammaCombMassFilter::hltFilter(edm::Event& iEvent,
0044                                         const edm::EventSetup& iSetup,
0045                                         trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0046   //right, issue 1, we dont know if this is a TriggerElectron, TriggerPhoton, TriggerCluster (should never be a TriggerCluster btw as that implies the 4-vectors are not stored in AOD)
0047 
0048   //trigger::TriggerObjectType firstLegTrigType;
0049   std::vector<math::XYZTLorentzVector> firstLegP4s;
0050 
0051   //trigger::TriggerObjectType secondLegTrigType;
0052   std::vector<math::XYZTLorentzVector> secondLegP4s;
0053 
0054   math::XYZTLorentzVector pairP4;
0055 
0056   getP4OfLegCands(iEvent, firstLegLastFilterToken_, firstLegP4s);
0057   getP4OfLegCands(iEvent, secondLegLastFilterToken_, secondLegP4s);
0058 
0059   bool accept = false;
0060   for (auto& firstLegP4 : firstLegP4s) {
0061     for (auto& secondLegP4 : secondLegP4s) {
0062       math::XYZTLorentzVector pairP4 = firstLegP4 + secondLegP4;
0063       double mass = pairP4.M();
0064       if (mass >= minMass_)
0065         accept = true;
0066     }
0067   }
0068 
0069   return accept;
0070 }
0071 
0072 void HLTEgammaCombMassFilter::getP4OfLegCands(const edm::Event& iEvent,
0073                                               const edm::EDGetTokenT<trigger::TriggerFilterObjectWithRefs>& filterToken,
0074                                               std::vector<math::XYZTLorentzVector>& p4s) {
0075   edm::Handle<trigger::TriggerFilterObjectWithRefs> filterOutput;
0076   iEvent.getByToken(filterToken, filterOutput);
0077 
0078   //its easier on the if statement flow if I try everything at once, shouldnt add to timing
0079   std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > phoCands;
0080   filterOutput->getObjects(trigger::TriggerPhoton, phoCands);
0081   std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > clusCands;
0082   filterOutput->getObjects(trigger::TriggerCluster, clusCands);
0083   std::vector<edm::Ref<reco::ElectronCollection> > eleCands;
0084   filterOutput->getObjects(trigger::TriggerElectron, eleCands);
0085 
0086   if (!phoCands.empty()) {  //its photons
0087     for (auto& phoCand : phoCands) {
0088       p4s.push_back(phoCand->p4());
0089     }
0090   } else if (!clusCands.empty()) {
0091     //try trigger cluster (should never be this, at the time of writing (17/1/11) this would indicate an error)
0092     for (auto& clusCand : clusCands) {
0093       p4s.push_back(clusCand->p4());
0094     }
0095   } else if (!eleCands.empty()) {
0096     for (auto& eleCand : eleCands) {
0097       p4s.push_back(eleCand->p4());
0098     }
0099   }
0100 }
0101 
0102 // declare this class as a framework plugin
0103 #include "FWCore/Framework/interface/MakerMacros.h"
0104 DEFINE_FWK_MODULE(HLTEgammaCombMassFilter);