Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-10-25 09:52:58

0001 /** \class HLTElectronGenericFilter
0002  *
0003  *
0004  *  \author Roberto Covarelli (CERN)
0005  *
0006  */
0007 
0008 #include "HLTElectronGenericFilter.h"
0009 
0010 #include "DataFormats/Common/interface/Handle.h"
0011 
0012 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0015 
0016 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
0017 #include "DataFormats/EgammaCandidates/interface/Electron.h"
0018 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
0019 #include "DataFormats/EgammaReco/interface/SuperClusterFwd.h"
0020 
0021 #include "DataFormats/Common/interface/AssociationMap.h"
0022 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
0023 
0024 //
0025 // constructors and destructor
0026 //
0027 HLTElectronGenericFilter::HLTElectronGenericFilter(const edm::ParameterSet& iConfig) : HLTFilter(iConfig) {
0028   candTag_ = iConfig.getParameter<edm::InputTag>("candTag");
0029   varTag_ = iConfig.getParameter<edm::InputTag>("varTag");
0030   lessThan_ = iConfig.getParameter<bool>("lessThan");
0031   thrRegularEB_ = iConfig.getParameter<double>("thrRegularEB");
0032   thrRegularEE_ = iConfig.getParameter<double>("thrRegularEE");
0033   thrOverPtEB_ = iConfig.getParameter<double>("thrOverPtEB");
0034   thrOverPtEE_ = iConfig.getParameter<double>("thrOverPtEE");
0035   thrTimesPtEB_ = iConfig.getParameter<double>("thrTimesPtEB");
0036   thrTimesPtEE_ = iConfig.getParameter<double>("thrTimesPtEE");
0037   ncandcut_ = iConfig.getParameter<int>("ncandcut");
0038   l1EGTag_ = iConfig.getParameter<edm::InputTag>("l1EGCand");
0039 
0040   candToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(candTag_);
0041   varToken_ = consumes<reco::ElectronIsolationMap>(varTag_);
0042 }
0043 
0044 HLTElectronGenericFilter::~HLTElectronGenericFilter() = default;
0045 
0046 void HLTElectronGenericFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0047   edm::ParameterSetDescription desc;
0048   makeHLTFilterDescription(desc);
0049   desc.add<edm::InputTag>("candTag", edm::InputTag("hltSingleElectronOneOEMinusOneOPFilter"));
0050   desc.add<edm::InputTag>("varTag", edm::InputTag("hltSingleElectronTrackIsol"));
0051   desc.add<bool>("lessThan", true);
0052   desc.add<double>("thrRegularEB", 0.0);
0053   desc.add<double>("thrRegularEE", 0.0);
0054   desc.add<double>("thrOverPtEB", -1.0);
0055   desc.add<double>("thrOverPtEE", -1.0);
0056   desc.add<double>("thrTimesPtEB", -1.0);
0057   desc.add<double>("thrTimesPtEE", -1.0);
0058   desc.add<int>("ncandcut", 1);
0059   desc.add<edm::InputTag>("l1EGCand", edm::InputTag("hltPixelMatchElectronsL1Iso"));
0060   descriptions.add("hltElectronGenericFilter", desc);
0061 }
0062 
0063 // ------------ method called to produce the data  ------------
0064 bool HLTElectronGenericFilter::hltFilter(edm::Event& iEvent,
0065                                          const edm::EventSetup& iSetup,
0066                                          trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0067   using namespace trigger;
0068 
0069   if (saveTags()) {
0070     filterproduct.addCollectionTag(l1EGTag_);
0071   }
0072 
0073   // Ref to Candidate object to be recorded in filter object
0074   reco::ElectronRef ref;
0075 
0076   edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
0077   iEvent.getByToken(candToken_, PrevFilterOutput);
0078 
0079   std::vector<edm::Ref<reco::ElectronCollection> > elecands;
0080   PrevFilterOutput->getObjects(TriggerElectron, elecands);
0081 
0082   //get hold of isolated association map
0083   edm::Handle<reco::ElectronIsolationMap> depMap;
0084   iEvent.getByToken(varToken_, depMap);
0085 
0086   // look at all photons, check cuts and add to filter object
0087   int n = 0;
0088 
0089   for (auto& elecand : elecands) {
0090     ref = elecand;
0091     reco::ElectronIsolationMap::const_iterator mapi = (*depMap).find(ref);
0092 
0093     float vali = mapi->val;
0094     float Pt = ref->pt();
0095     float Eta = fabs(ref->eta());
0096 
0097     if (lessThan_) {
0098       if ((Eta < 1.479 && vali <= thrRegularEB_) || (Eta >= 1.479 && vali <= thrRegularEE_)) {
0099         n++;
0100         filterproduct.addObject(TriggerElectron, ref);
0101         continue;
0102       }
0103       if (Pt > 0. && (thrOverPtEB_ > 0. || thrOverPtEE_ > 0. || thrTimesPtEB_ > 0. || thrTimesPtEE_ > 0.)) {
0104         if ((Eta < 1.479 && vali / Pt <= thrOverPtEB_) || (Eta >= 1.479 && vali / Pt <= thrOverPtEE_)) {
0105           n++;
0106           filterproduct.addObject(TriggerElectron, ref);
0107           continue;
0108         }
0109         if ((Eta < 1.479 && vali * Pt <= thrTimesPtEB_) || (Eta >= 1.479 && vali * Pt <= thrTimesPtEE_)) {
0110           n++;
0111           filterproduct.addObject(TriggerElectron, ref);
0112         }
0113       }
0114     } else {
0115       if ((Eta < 1.479 && vali >= thrRegularEB_) || (Eta >= 1.479 && vali >= thrRegularEE_)) {
0116         n++;
0117         filterproduct.addObject(TriggerElectron, ref);
0118         continue;
0119       }
0120       if (Pt > 0. && (thrOverPtEB_ > 0. || thrOverPtEE_ > 0. || thrTimesPtEB_ > 0. || thrTimesPtEE_ > 0.)) {
0121         if ((Eta < 1.479 && vali / Pt >= thrOverPtEB_) || (Eta >= 1.479 && vali / Pt >= thrOverPtEE_)) {
0122           n++;
0123           filterproduct.addObject(TriggerElectron, ref);
0124           continue;
0125         }
0126         if ((Eta < 1.479 && vali * Pt >= thrTimesPtEB_) || (Eta >= 1.479 && vali * Pt >= thrTimesPtEE_)) {
0127           n++;
0128           filterproduct.addObject(TriggerElectron, ref);
0129         }
0130       }
0131     }
0132   }
0133 
0134   // filter decision
0135   bool accept(n >= ncandcut_);
0136 
0137   return accept;
0138 }
0139 
0140 // declare this class as a framework plugin
0141 #include "FWCore/Framework/interface/MakerMacros.h"
0142 DEFINE_FWK_MODULE(HLTElectronGenericFilter);