Back to home page

Project CMSSW displayed by LXR

 
 

    


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

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