Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:27:44

0001 // -*- C++ -*-
0002 //
0003 // Package:     RecoTauTag/HLTProducers
0004 // Class  :     TrackingRegionsFromBeamSpotAndL2Tau
0005 //
0006 // Implementation:
0007 //     [Notes on implementation]
0008 //
0009 // Original Author:  Christopher Jones
0010 //         Created:  Tue, 13 Sep 2022 21:13:17 GMT
0011 //
0012 
0013 // system include files
0014 
0015 // user include files
0016 #include "TrackingRegionsFromBeamSpotAndL2Tau.h"
0017 
0018 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegionProducerFactory.h"
0019 
0020 DEFINE_EDM_PLUGIN(TrackingRegionProducerFactory,
0021                   TrackingRegionsFromBeamSpotAndL2Tau,
0022                   "TrackingRegionsFromBeamSpotAndL2Tau");
0023 
0024 TrackingRegionsFromBeamSpotAndL2Tau::TrackingRegionsFromBeamSpotAndL2Tau(const edm::ParameterSet& conf,
0025                                                                          edm::ConsumesCollector&& iC) {
0026   edm::LogInfo("TrackingRegionsFromBeamSpotAndL2Tau") << "Enter the TrackingRegionsFromBeamSpotAndL2Tau";
0027 
0028   edm::ParameterSet regionPSet = conf.getParameter<edm::ParameterSet>("RegionPSet");
0029 
0030   m_ptMin = regionPSet.getParameter<double>("ptMin");
0031   m_originRadius = regionPSet.getParameter<double>("originRadius");
0032   m_originHalfLength = regionPSet.getParameter<double>("originHalfLength");
0033   m_deltaEta = regionPSet.getParameter<double>("deltaEta");
0034   m_deltaPhi = regionPSet.getParameter<double>("deltaPhi");
0035   token_jet = iC.consumes<reco::CandidateView>(regionPSet.getParameter<edm::InputTag>("JetSrc"));
0036   m_jetMinPt = regionPSet.getParameter<double>("JetMinPt");
0037   m_jetMaxEta = regionPSet.getParameter<double>("JetMaxEta");
0038   m_jetMaxN = regionPSet.getParameter<int>("JetMaxN");
0039   token_beamSpot = iC.consumes<reco::BeamSpot>(regionPSet.getParameter<edm::InputTag>("beamSpot"));
0040   m_precise = regionPSet.getParameter<bool>("precise");
0041 
0042   if (regionPSet.exists("searchOpt"))
0043     m_searchOpt = regionPSet.getParameter<bool>("searchOpt");
0044   else
0045     m_searchOpt = false;
0046 
0047   m_whereToUseMeasurementTracker = RectangularEtaPhiTrackingRegion::stringToUseMeasurementTracker(
0048       regionPSet.getParameter<std::string>("whereToUseMeasurementTracker"));
0049   if (m_whereToUseMeasurementTracker != RectangularEtaPhiTrackingRegion::UseMeasurementTracker::kNever) {
0050     token_measurementTracker =
0051         iC.consumes<MeasurementTrackerEvent>(regionPSet.getParameter<edm::InputTag>("measurementTrackerName"));
0052   }
0053   token_field = iC.esConsumes();
0054   if (m_precise) {
0055     token_msmaker = iC.esConsumes();
0056   }
0057 }
0058 
0059 void TrackingRegionsFromBeamSpotAndL2Tau::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0060   edm::ParameterSetDescription desc;
0061 
0062   desc.add<double>("ptMin", 5.0);
0063   desc.add<double>("originRadius", 0.2);
0064   desc.add<double>("originHalfLength", 24.0);
0065   desc.add<double>("deltaEta", 0.3);
0066   desc.add<double>("deltaPhi", 0.3);
0067   desc.add<edm::InputTag>("JetSrc", edm::InputTag("hltFilterL2EtCutDoublePFIsoTau25Trk5"));
0068   desc.add<double>("JetMinPt", 25.0);
0069   desc.add<double>("JetMaxEta", 2.1);
0070   desc.add<int>("JetMaxN", 10);
0071   desc.add<edm::InputTag>("beamSpot", edm::InputTag("hltOnlineBeamSpot"));
0072   desc.add<bool>("precise", true);
0073   desc.add<std::string>("howToUseMeasurementTracker", "Never");
0074   desc.add<edm::InputTag>("measurementTrackerName", edm::InputTag("MeasurementTrackerEvent"));
0075 
0076   // Only for backwards-compatibility
0077   edm::ParameterSetDescription descRegion;
0078   descRegion.add<edm::ParameterSetDescription>("RegionPSet", desc);
0079 
0080   descriptions.add("trackingRegionsFromBeamSpotAndL2Tau", descRegion);
0081 }
0082 
0083 std::vector<std::unique_ptr<TrackingRegion> > TrackingRegionsFromBeamSpotAndL2Tau::regions(
0084     const edm::Event& e, const edm::EventSetup& es) const {
0085   std::vector<std::unique_ptr<TrackingRegion> > result;
0086 
0087   // use beam spot to pick up the origin
0088   edm::Handle<reco::BeamSpot> bsHandle;
0089   e.getByToken(token_beamSpot, bsHandle);
0090   if (!bsHandle.isValid())
0091     return result;
0092   const reco::BeamSpot& bs = *bsHandle;
0093   GlobalPoint origin(bs.x0(), bs.y0(), bs.z0());
0094 
0095   // pick up the candidate objects of interest
0096   edm::Handle<reco::CandidateView> objects;
0097   e.getByToken(token_jet, objects);
0098   size_t n_objects = objects->size();
0099   if (n_objects == 0)
0100     return result;
0101 
0102   const MeasurementTrackerEvent* measurementTracker = nullptr;
0103   if (!token_measurementTracker.isUninitialized()) {
0104     edm::Handle<MeasurementTrackerEvent> hmte;
0105     e.getByToken(token_measurementTracker, hmte);
0106     measurementTracker = hmte.product();
0107   }
0108 
0109   const auto& field = es.getData(token_field);
0110   const MultipleScatteringParametrisationMaker* msmaker = nullptr;
0111   if (m_precise) {
0112     msmaker = &es.getData(token_msmaker);
0113   }
0114 
0115   // create maximum JetMaxN tracking regions in directions of
0116   // highest pt jets that are above threshold and are within allowed eta
0117   // (we expect that jet collection was sorted in decreasing pt order)
0118   int n_regions = 0;
0119   for (size_t i = 0; i < n_objects && n_regions < m_jetMaxN; ++i) {
0120     const reco::Candidate& jet = (*objects)[i];
0121     if (jet.pt() < m_jetMinPt || std::abs(jet.eta()) > m_jetMaxEta)
0122       continue;
0123 
0124     GlobalVector direction(jet.momentum().x(), jet.momentum().y(), jet.momentum().z());
0125 
0126     result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(direction,
0127                                                                        origin,
0128                                                                        m_ptMin,
0129                                                                        m_originRadius,
0130                                                                        m_originHalfLength,
0131                                                                        m_deltaEta,
0132                                                                        m_deltaPhi,
0133                                                                        field,
0134                                                                        msmaker,
0135                                                                        m_precise,
0136                                                                        m_whereToUseMeasurementTracker,
0137                                                                        measurementTracker,
0138                                                                        m_searchOpt));
0139     ++n_regions;
0140   }
0141   //std::cout<<"nregions = "<<n_regions<<std::endl;
0142   return result;
0143 }