CMS 3D CMS Logo

TauRegionalPixelSeedGenerator.h
Go to the documentation of this file.
1 #ifndef TauRegionalPixelSeedGenerator_h
2 #define TauRegionalPixelSeedGenerator_h
3 
4 //
5 // Class: TauRegionalPixelSeedGenerator
6 
27 // Math
28 #include "Math/GenVector/VectorUtil.h"
29 #include "Math/GenVector/PxPyPzE4D.h"
30 
32 
35 
37 public:
40  edm::LogInfo("TauRegionalPixelSeedGenerator") << "Enter the TauRegionalPixelSeedGenerator";
41 
42  edm::ParameterSet regionPSet = conf.getParameter<edm::ParameterSet>("RegionPSet");
43 
44  m_ptMin = regionPSet.getParameter<double>("ptMin");
45  m_originRadius = regionPSet.getParameter<double>("originRadius");
46  m_halfLength = regionPSet.getParameter<double>("originHalfLength");
47  m_deltaEta = regionPSet.getParameter<double>("deltaEtaRegion");
48  m_deltaPhi = regionPSet.getParameter<double>("deltaPhiRegion");
49  token_jet = iC.consumes<reco::CandidateView>(regionPSet.getParameter<edm::InputTag>("JetSrc"));
50  token_vertex = iC.consumes<reco::VertexCollection>(regionPSet.getParameter<edm::InputTag>("vertexSrc"));
51  if (regionPSet.exists("searchOpt")) {
52  m_searchOpt = regionPSet.getParameter<bool>("searchOpt");
53  } else {
54  m_searchOpt = false;
55  }
57  regionPSet.getParameter<std::string>("howToUseMeasurementTracker"));
59  // temporary until everything migrated to InputTag
61  if (regionPSet.existsAs<edm::InputTag>("measurementTrackerName")) {
62  tag = regionPSet.getParameter<edm::InputTag>("measurementTrackerName");
63  } else {
64  tag = edm::InputTag(regionPSet.getParameter<std::string>("measurementTrackerName"));
65  }
66 
68  }
69  }
70 
72 
73  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
75 
76  desc.add<double>("ptMin", 5.0);
77  desc.add<double>("originHalfLength", 0.2);
78  desc.add<double>("originRadius", 0.2);
79  desc.add<double>("deltaEtaRegion", 0.1);
80  desc.add<double>("deltaPhiRegion", 0.1);
81  desc.add<edm::InputTag>("JetSrc", edm::InputTag("icone5Tau1"));
82  desc.add<edm::InputTag>("vertexSrc", edm::InputTag("pixelVertices"));
83  desc.add<bool>("searchOpt", false);
84 
85  desc.add<std::string>("howToUseMeasurementTracker", "ForSiStrips");
86 
87  // allow both InputTag and string for the moment, use InputTag as the default
89  "measurementTrackerName", edm::InputTag("MeasurementTrackerEvent"), true) xor
90  edm::ParameterDescription<std::string>("measurementTrackerName", "MeasurementTrackerEvent", true));
91 
92  // Only for backwards-compatibility
94  descRegion.add<edm::ParameterSetDescription>("RegionPSet", desc);
95 
96  descriptions.add("tauRegionalPixelSeedTrackingRegions", descRegion);
97  }
98 
99  std::vector<std::unique_ptr<TrackingRegion> > regions(const edm::Event& e, const edm::EventSetup& es) const override {
100  std::vector<std::unique_ptr<TrackingRegion> > result;
101 
102  // double originZ;
103  double deltaZVertex, deltaRho;
105  // get the primary vertex
107  e.getByToken(token_vertex, h_vertices);
108  const reco::VertexCollection& vertices = *h_vertices;
109  if (not vertices.empty()) {
110  // originZ = vertices.front().z();
111  GlobalPoint myTmp(vertices.at(0).position().x(), vertices.at(0).position().y(), vertices.at(0).position().z());
112  vertex = myTmp;
113  deltaZVertex = m_halfLength;
114  deltaRho = m_originRadius;
115  } else {
116  // originZ = 0.;
117  GlobalPoint myTmp(0., 0., 0.);
118  vertex = myTmp;
119  deltaRho = 1.;
120  deltaZVertex = 15.;
121  }
122 
123  // get the jet direction
125  e.getByToken(token_jet, h_jets);
126 
130  e.getByToken(token_measurementTracker, hmte);
131  measurementTracker = hmte.product();
132  }
133 
134  const auto& bfield = es.getData(token_bfield);
135  const auto& msmaker = es.getData(token_msmaker);
136 
137  for (const reco::Candidate& myJet : *h_jets) {
138  GlobalVector jetVector(myJet.momentum().x(), myJet.momentum().y(), myJet.momentum().z());
139  // GlobalPoint vertex(0, 0, originZ);
140  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(jetVector,
141  vertex,
142  m_ptMin,
143  deltaRho,
144  deltaZVertex,
145  m_deltaEta,
146  m_deltaPhi,
147  bfield,
148  &msmaker,
149  true,
152  m_searchOpt));
153  }
154 
155  return result;
156  }
157 
158 private:
159  float m_ptMin;
162  float m_deltaEta;
163  float m_deltaPhi;
171 };
172 
173 #endif
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
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:104
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
static UseMeasurementTracker stringToUseMeasurementTracker(const std::string &name)
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:171
edm::EDGetTokenT< MeasurementTrackerEvent > token_measurementTracker
TauRegionalPixelSeedGenerator(const edm::ParameterSet &conf, edm::ConsumesCollector &&iC)
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > token_bfield
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< reco::VertexCollection > token_vertex
Log< level::Info, false > LogInfo
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_howToUseMeasurementTracker
std::vector< std::unique_ptr< TrackingRegion > > regions(const edm::Event &e, const edm::EventSetup &es) const override
edm::ESGetToken< MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord > token_msmaker
edm::EDGetTokenT< reco::CandidateView > token_jet