Back to home page

Project CMSSW displayed by LXR

 
 

    


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

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