CMS 3D CMS Logo

TrackingRegionsFromBeamSpotAndL2Tau.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: RecoTauTag/HLTProducers
4 // Class : TrackingRegionsFromBeamSpotAndL2Tau
5 //
6 // Implementation:
7 // [Notes on implementation]
8 //
9 // Original Author: Christopher Jones
10 // Created: Tue, 13 Sep 2022 21:13:17 GMT
11 //
12 
13 // system include files
14 
15 // user include files
17 
19 
22  "TrackingRegionsFromBeamSpotAndL2Tau");
23 
26  edm::LogInfo("TrackingRegionsFromBeamSpotAndL2Tau") << "Enter the TrackingRegionsFromBeamSpotAndL2Tau";
27 
28  edm::ParameterSet regionPSet = conf.getParameter<edm::ParameterSet>("RegionPSet");
29 
30  m_ptMin = regionPSet.getParameter<double>("ptMin");
31  m_originRadius = regionPSet.getParameter<double>("originRadius");
32  m_originHalfLength = regionPSet.getParameter<double>("originHalfLength");
33  m_deltaEta = regionPSet.getParameter<double>("deltaEta");
34  m_deltaPhi = regionPSet.getParameter<double>("deltaPhi");
35  token_jet = iC.consumes<reco::CandidateView>(regionPSet.getParameter<edm::InputTag>("JetSrc"));
36  m_jetMinPt = regionPSet.getParameter<double>("JetMinPt");
37  m_jetMaxEta = regionPSet.getParameter<double>("JetMaxEta");
38  m_jetMaxN = regionPSet.getParameter<int>("JetMaxN");
39  token_beamSpot = iC.consumes<reco::BeamSpot>(regionPSet.getParameter<edm::InputTag>("beamSpot"));
40  m_precise = regionPSet.getParameter<bool>("precise");
41 
42  if (regionPSet.exists("searchOpt"))
43  m_searchOpt = regionPSet.getParameter<bool>("searchOpt");
44  else
45  m_searchOpt = false;
46 
48  regionPSet.getParameter<std::string>("whereToUseMeasurementTracker"));
51  iC.consumes<MeasurementTrackerEvent>(regionPSet.getParameter<edm::InputTag>("measurementTrackerName"));
52  }
53  token_field = iC.esConsumes();
54  if (m_precise) {
55  token_msmaker = iC.esConsumes();
56  }
57 }
58 
61 
62  desc.add<double>("ptMin", 5.0);
63  desc.add<double>("originRadius", 0.2);
64  desc.add<double>("originHalfLength", 24.0);
65  desc.add<double>("deltaEta", 0.3);
66  desc.add<double>("deltaPhi", 0.3);
67  desc.add<edm::InputTag>("JetSrc", edm::InputTag("hltFilterL2EtCutDoublePFIsoTau25Trk5"));
68  desc.add<double>("JetMinPt", 25.0);
69  desc.add<double>("JetMaxEta", 2.1);
70  desc.add<int>("JetMaxN", 10);
71  desc.add<edm::InputTag>("beamSpot", edm::InputTag("hltOnlineBeamSpot"));
72  desc.add<bool>("precise", true);
73  desc.add<std::string>("howToUseMeasurementTracker", "Never");
74  desc.add<edm::InputTag>("measurementTrackerName", edm::InputTag("MeasurementTrackerEvent"));
75 
76  // Only for backwards-compatibility
78  descRegion.add<edm::ParameterSetDescription>("RegionPSet", desc);
79 
80  descriptions.add("trackingRegionsFromBeamSpotAndL2Tau", descRegion);
81 }
82 
83 std::vector<std::unique_ptr<TrackingRegion> > TrackingRegionsFromBeamSpotAndL2Tau::regions(
84  const edm::Event& e, const edm::EventSetup& es) const {
85  std::vector<std::unique_ptr<TrackingRegion> > result;
86 
87  // use beam spot to pick up the origin
89  e.getByToken(token_beamSpot, bsHandle);
90  if (!bsHandle.isValid())
91  return result;
92  const reco::BeamSpot& bs = *bsHandle;
93  GlobalPoint origin(bs.x0(), bs.y0(), bs.z0());
94 
95  // pick up the candidate objects of interest
97  e.getByToken(token_jet, objects);
98  size_t n_objects = objects->size();
99  if (n_objects == 0)
100  return result;
101 
105  e.getByToken(token_measurementTracker, hmte);
106  measurementTracker = hmte.product();
107  }
108 
109  const auto& field = es.getData(token_field);
110  const MultipleScatteringParametrisationMaker* msmaker = nullptr;
111  if (m_precise) {
112  msmaker = &es.getData(token_msmaker);
113  }
114 
115  // create maximum JetMaxN tracking regions in directions of
116  // highest pt jets that are above threshold and are within allowed eta
117  // (we expect that jet collection was sorted in decreasing pt order)
118  int n_regions = 0;
119  for (size_t i = 0; i < n_objects && n_regions < m_jetMaxN; ++i) {
120  const reco::Candidate& jet = (*objects)[i];
121  if (jet.pt() < m_jetMinPt || std::abs(jet.eta()) > m_jetMaxEta)
122  continue;
123 
124  GlobalVector direction(jet.momentum().x(), jet.momentum().y(), jet.momentum().z());
125 
126  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(direction,
127  origin,
128  m_ptMin,
131  m_deltaEta,
132  m_deltaPhi,
133  field,
134  msmaker,
135  m_precise,
138  m_searchOpt));
139  ++n_regions;
140  }
141  //std::cout<<"nregions = "<<n_regions<<std::endl;
142  return result;
143 }
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_whereToUseMeasurementTracker
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
edm::ESGetToken< MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord > token_msmaker
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool exists(std::string const &parameterName) const
checks if a parameter exists
T const * product() const
Definition: Handle.h:70
constexpr bool isUninitialized() const noexcept
Definition: EDGetToken.h:98
edm::EDGetTokenT< MeasurementTrackerEvent > token_measurementTracker
static UseMeasurementTracker stringToUseMeasurementTracker(const std::string &name)
edm::EDGetTokenT< reco::CandidateView > token_jet
std::vector< std::unique_ptr< TrackingRegion > > regions(const edm::Event &e, const edm::EventSetup &es) const override
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > token_field
Log< level::Info, false > LogInfo
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool isValid() const
Definition: HandleBase.h:70
#define DEFINE_EDM_PLUGIN(factory, type, name)
TrackingRegionsFromBeamSpotAndL2Tau(const edm::ParameterSet &conf, edm::ConsumesCollector &&iC)