Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:18:19

0001 /** \class HLTDisplacedEgammaFilter
0002  *
0003  *
0004  *  \author Monica Vazquez Acosta (CERN)
0005  *
0006  */
0007 #include <cmath>
0008 
0009 #include "HLTDisplacedEgammaFilter.h"
0010 #include "DataFormats/Common/interface/Handle.h"
0011 #include "DataFormats/Math/interface/deltaR.h"
0012 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
0013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0014 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0015 #include "RecoEcal/EgammaCoreTools/interface/EcalClusterTools.h"
0016 
0017 HLTDisplacedEgammaFilter::HLTDisplacedEgammaFilter(const edm::ParameterSet& iConfig) : HLTFilter(iConfig) {
0018   inputTag_ = iConfig.getParameter<edm::InputTag>("inputTag");
0019   ncandcut_ = iConfig.getParameter<int>("ncandcut");
0020   l1EGTag_ = iConfig.getParameter<edm::InputTag>("l1EGCand");
0021 
0022   inputTrk = iConfig.getParameter<edm::InputTag>("inputTrack");
0023   trkPtCut = iConfig.getParameter<double>("trackPtCut");
0024 
0025   // track dR^2 threshold with sign
0026   auto const trkDrCut = iConfig.getParameter<double>("trackdRCut");
0027   trkDr2Cut = trkDrCut * std::abs(trkDrCut);
0028 
0029   maxTrkCut = iConfig.getParameter<int>("maxTrackCut");
0030 
0031   rechitsEB = iConfig.getParameter<edm::InputTag>("RecHitsEB");
0032   rechitsEE = iConfig.getParameter<edm::InputTag>("RecHitsEE");
0033 
0034   EBOnly = iConfig.getParameter<bool>("EBOnly");
0035   sMin_min = iConfig.getParameter<double>("sMin_min");
0036   sMin_max = iConfig.getParameter<double>("sMin_max");
0037   sMaj_min = iConfig.getParameter<double>("sMaj_min");
0038   sMaj_max = iConfig.getParameter<double>("sMaj_max");
0039   seedTimeMin = iConfig.getParameter<double>("seedTimeMin");
0040   seedTimeMax = iConfig.getParameter<double>("seedTimeMax");
0041   useTrackVeto = iConfig.getParameter<bool>("useTrackVeto");
0042 
0043   inputToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(inputTag_);
0044   rechitsEBToken_ = consumes<EcalRecHitCollection>(rechitsEB);
0045   rechitsEEToken_ = consumes<EcalRecHitCollection>(rechitsEE);
0046 
0047   if (useTrackVeto) {
0048     inputTrkToken_ = consumes<reco::TrackCollection>(inputTrk);
0049   }
0050 }
0051 
0052 void HLTDisplacedEgammaFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0053   edm::ParameterSetDescription desc;
0054   makeHLTFilterDescription(desc);
0055   desc.add<edm::InputTag>("inputTag", edm::InputTag("hltEGRegionalL1SingleEG22"));
0056   desc.add<edm::InputTag>("l1EGCand", edm::InputTag("hltL1IsoRecoEcalCandidate"));
0057   desc.add<edm::InputTag>("RecHitsEB", edm::InputTag("hltEcalRecHit", "EcalRecHitsEB"));
0058   desc.add<edm::InputTag>("RecHitsEE", edm::InputTag("hltEcalRecHit", "EcalRecHitsEE"));
0059   desc.add<edm::InputTag>("inputTrack", edm::InputTag("hltL1SeededEgammaRegionalCTFFinalFitWithMaterial"));
0060   desc.add<int>("ncandcut", 1);
0061   desc.add<bool>("EBOnly", false);
0062   desc.add<double>("sMin_min", 0.1);
0063   desc.add<double>("sMin_max", 0.4);
0064   desc.add<double>("sMaj_min", 0.0);
0065   desc.add<double>("sMaj_max", 999.0);
0066   desc.add<double>("seedTimeMin", -25.0);
0067   desc.add<double>("seedTimeMax", 25.0);
0068   desc.add<bool>("useTrackVeto", true);
0069   desc.add<int>("maxTrackCut", 0);
0070   desc.add<double>("trackPtCut", 3.0);
0071   desc.add<double>("trackdRCut", 0.5);
0072   descriptions.add("hltDisplacedEgammaFilter", desc);
0073 }
0074 
0075 // ------------ method called to produce the data  ------------
0076 
0077 bool HLTDisplacedEgammaFilter::hltFilter(edm::Event& iEvent,
0078                                          const edm::EventSetup& iSetup,
0079                                          trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0080   using namespace trigger;
0081 
0082   // The filter object
0083   if (saveTags()) {
0084     filterproduct.addCollectionTag(l1EGTag_);
0085   }
0086 
0087   // get hold of filtered candidates
0088   //edm::Handle<reco::HLTFilterObjectWithRefs> recoecalcands;
0089   edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
0090   iEvent.getByToken(inputToken_, PrevFilterOutput);
0091 
0092   // get hold of collection of objects
0093   edm::Handle<reco::TrackCollection> tracks;
0094   if (useTrackVeto) {
0095     iEvent.getByToken(inputTrkToken_, tracks);
0096   }
0097 
0098   // get the EcalRecHit
0099   edm::Handle<EcalRecHitCollection> rechitsEB_;
0100   edm::Handle<EcalRecHitCollection> rechitsEE_;
0101   iEvent.getByToken(rechitsEBToken_, rechitsEB_);
0102   iEvent.getByToken(rechitsEEToken_, rechitsEE_);
0103 
0104   std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > recoecalcands;
0105   PrevFilterOutput->getObjects(TriggerCluster, recoecalcands);
0106   if (recoecalcands.empty())
0107     PrevFilterOutput->getObjects(TriggerPhoton, recoecalcands);
0108 
0109   // look at all candidates,  check cuts and add to filter object
0110   int n(0);
0111 
0112   for (auto const& ref : recoecalcands) {
0113     if (EBOnly && std::abs(ref->eta()) >= 1.479)
0114       continue;
0115 
0116     // S_Minor Cuts from the seed cluster
0117     reco::CaloClusterPtr SCseed = ref->superCluster()->seed();
0118     const EcalRecHitCollection* rechits = (std::abs(ref->eta()) < 1.479) ? rechitsEB_.product() : rechitsEE_.product();
0119 
0120     Cluster2ndMoments moments = EcalClusterTools::cluster2ndMoments(*SCseed, *rechits);
0121     float sMin = moments.sMin;
0122     float sMaj = moments.sMaj;
0123     if (sMin < sMin_min || sMin > sMin_max)
0124       continue;
0125     if (sMaj < sMaj_min || sMaj > sMaj_max)
0126       continue;
0127 
0128     // seed Time
0129     std::pair<DetId, float> maxRH = EcalClusterTools::getMaximum(*SCseed, rechits);
0130     DetId seedCrystalId = maxRH.first;
0131     auto seedRH = rechits->find(seedCrystalId);
0132     float seedTime = (float)seedRH->time();
0133     if (seedTime < seedTimeMin || seedTime > seedTimeMax)
0134       continue;
0135 
0136     //Track Veto
0137 
0138     if (useTrackVeto) {
0139       int nTrk = 0;
0140       for (auto const& it : *tracks) {
0141         if (it.pt() < trkPtCut)
0142           continue;
0143         auto const dR2 = reco::deltaR2(it.eta(), it.phi(), ref->eta(), ref->phi());
0144         if (dR2 < trkDr2Cut)
0145           nTrk++;
0146         if (nTrk > maxTrkCut)
0147           break;
0148       }
0149       if (nTrk > maxTrkCut)
0150         continue;
0151     }
0152 
0153     n++;
0154     // std::cout << "Passed eta: " << ref->eta() << std::endl;
0155     filterproduct.addObject(TriggerCluster, ref);
0156   }
0157 
0158   // filter decision
0159   return (n >= ncandcut_);
0160 }
0161 
0162 // declare this class as a framework plugin
0163 #include "FWCore/Framework/interface/MakerMacros.h"
0164 DEFINE_FWK_MODULE(HLTDisplacedEgammaFilter);