Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-07 04:36:36

0001 /** \class HLTElectronEoverpFilterRegional
0002  *
0003  *  \author Monica Vazquez Acosta (CERN)
0004  *
0005  *
0006  */
0007 
0008 #include "HLTElectronEoverpFilterRegional.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/RecoCandidate/interface/RecoEcalCandidateFwd.h"
0018 #include "DataFormats/TrackReco/interface/Track.h"
0019 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
0020 
0021 //
0022 // constructors and destructor
0023 //
0024 HLTElectronEoverpFilterRegional::HLTElectronEoverpFilterRegional(const edm::ParameterSet& iConfig)
0025     : HLTFilter(iConfig) {
0026   candTag_ = iConfig.getParameter<edm::InputTag>("candTag");
0027   electronIsolatedProducer_ = iConfig.getParameter<edm::InputTag>("electronIsolatedProducer");
0028   electronNonIsolatedProducer_ = iConfig.getParameter<edm::InputTag>("electronNonIsolatedProducer");
0029   eoverpbarrelcut_ = iConfig.getParameter<double>("eoverpbarrelcut");
0030   eoverpendcapcut_ = iConfig.getParameter<double>("eoverpendcapcut");
0031   ncandcut_ = iConfig.getParameter<int>("ncandcut");
0032   doIsolated_ = iConfig.getParameter<bool>("doIsolated");
0033   candToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(candTag_);
0034   electronIsolatedToken_ = consumes<reco::ElectronCollection>(electronIsolatedProducer_);
0035   if (!doIsolated_)
0036     electronNonIsolatedToken_ = consumes<reco::ElectronCollection>(electronNonIsolatedProducer_);
0037 }
0038 
0039 HLTElectronEoverpFilterRegional::~HLTElectronEoverpFilterRegional() = default;
0040 
0041 void HLTElectronEoverpFilterRegional::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0042   edm::ParameterSetDescription desc;
0043   makeHLTFilterDescription(desc);
0044   desc.add<edm::InputTag>("candTag", edm::InputTag("hltElectronPixelMatchFilter"));
0045   desc.add<edm::InputTag>("electronIsolatedProducer", edm::InputTag("pixelMatchElectronsForHLT"));
0046   desc.add<edm::InputTag>("electronNonIsolatedProducer", edm::InputTag("pixelMatchElectronsForHLT"));
0047   desc.add<double>("eoverpbarrelcut", 1.5);
0048   desc.add<double>("eoverpendcapcut", 2.45);
0049   desc.add<int>("ncandcut", 1);
0050   desc.add<bool>("doIsolated", true);
0051   descriptions.add("hltElectronEoverpFilter", desc);
0052 }
0053 
0054 // ------------ method called to produce the data  ------------
0055 bool HLTElectronEoverpFilterRegional::hltFilter(edm::Event& iEvent,
0056                                                 const edm::EventSetup& iSetup,
0057                                                 trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0058   // The filter object
0059   using namespace trigger;
0060   if (saveTags()) {
0061     filterproduct.addCollectionTag(electronIsolatedProducer_);
0062     if (not doIsolated_)
0063       filterproduct.addCollectionTag(electronNonIsolatedProducer_);
0064   }
0065 
0066   edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
0067   iEvent.getByToken(candToken_, PrevFilterOutput);
0068 
0069   std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > recoecalcands;
0070   PrevFilterOutput->getObjects(TriggerCluster, recoecalcands);
0071   if (recoecalcands.empty())
0072     PrevFilterOutput->getObjects(TriggerPhoton, recoecalcands);
0073   // Get the HLT electrons from EgammaHLTPixelMatchElectronProducers
0074   edm::Handle<reco::ElectronCollection> electronIsolatedHandle;
0075   iEvent.getByToken(electronIsolatedToken_, electronIsolatedHandle);
0076 
0077   edm::Handle<reco::ElectronCollection> electronNonIsolatedHandle;
0078   if (!doIsolated_) {
0079     iEvent.getByToken(electronNonIsolatedToken_, electronNonIsolatedHandle);
0080   }
0081 
0082   // look at all candidates,  check cuts and add to filter object
0083   int n(0);
0084 
0085   //loop over all the RecoCandidates from the previous filter,
0086   // associate them with the corresponding Electron object
0087   //(the matching is done checking the super clusters)
0088   // and put into the event a Ref to the Electron objects that passes the
0089   // selections
0090   for (auto& recoecalcand : recoecalcands) {
0091     reco::SuperClusterRef recr2 = recoecalcand->superCluster();
0092 
0093     //loop over the electrons to find the matching one
0094     for (auto iElectron = electronIsolatedHandle->begin(); iElectron != electronIsolatedHandle->end(); iElectron++) {
0095       reco::ElectronRef electronref(
0096           reco::ElectronRef(electronIsolatedHandle, iElectron - electronIsolatedHandle->begin()));
0097       const reco::SuperClusterRef theClus = electronref->superCluster();
0098 
0099       if (&(*recr2) == &(*theClus)) {
0100         float elecEoverp = 0;
0101         const math::XYZVector trackMom = electronref->track()->momentum();
0102         if (trackMom.R() != 0)
0103           elecEoverp = electronref->superCluster()->energy() / trackMom.R();
0104 
0105         if (fabs(electronref->eta()) < 1.5) {
0106           if (elecEoverp < eoverpbarrelcut_) {
0107             n++;
0108             filterproduct.addObject(TriggerElectron, electronref);
0109           }
0110         }
0111         if (fabs(electronref->eta()) > 1.5) {
0112           if (elecEoverp < eoverpendcapcut_) {
0113             n++;
0114             filterproduct.addObject(TriggerElectron, electronref);
0115           }
0116         }
0117       }  //end of the if checking the matching of the SC from RecoCandidate and the one from Electrons
0118     }  //end of loop over electrons
0119 
0120     if (!doIsolated_) {
0121       //loop over the electrons to find the matching one
0122       for (auto iElectron = electronNonIsolatedHandle->begin(); iElectron != electronNonIsolatedHandle->end();
0123            iElectron++) {
0124         reco::ElectronRef electronref(
0125             reco::ElectronRef(electronNonIsolatedHandle, iElectron - electronNonIsolatedHandle->begin()));
0126         const reco::SuperClusterRef theClus = electronref->superCluster();
0127 
0128         if (&(*recr2) == &(*theClus)) {
0129           float elecEoverp = 0;
0130           const math::XYZVector trackMom = electronref->track()->momentum();
0131           if (trackMom.R() != 0)
0132             elecEoverp = electronref->superCluster()->energy() / trackMom.R();
0133 
0134           if (fabs(electronref->eta()) < 1.5) {
0135             if (elecEoverp < eoverpbarrelcut_) {
0136               n++;
0137               filterproduct.addObject(TriggerElectron, electronref);
0138             }
0139           }
0140           if (fabs(electronref->eta()) > 1.5) {
0141             if (elecEoverp < eoverpendcapcut_) {
0142               n++;
0143               filterproduct.addObject(TriggerElectron, electronref);
0144             }
0145           }
0146         }  //end of the if checking the matching of the SC from RecoCandidate and the one from Electrons
0147       }  //end of loop over electrons
0148     }
0149   }  //end of loop ober candidates
0150 
0151   // filter decision
0152   bool accept(n >= ncandcut_);
0153 
0154   return accept;
0155 }
0156 
0157 // declare this class as a framework plugin
0158 #include "FWCore/Framework/interface/MakerMacros.h"
0159 DEFINE_FWK_MODULE(HLTElectronEoverpFilterRegional);