CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HITRegionalPixelSeedGenerator.cc
Go to the documentation of this file.
1 //
2 // Class: HITRegionalPixelSeedGenerator
3 
24 
25 // Math
26 #include "Math/GenVector/VectorUtil.h"
27 #include "Math/GenVector/PxPyPzE4D.h"
28 
30 
33 
35 public:
37  : m_regionPSet(conf.getParameter<edm::ParameterSet>("RegionPSet")),
38  ptmin(m_regionPSet.getParameter<double>("ptMin")),
39  originradius(m_regionPSet.getParameter<double>("originRadius")),
40  halflength(m_regionPSet.getParameter<double>("originHalfLength")),
41  etaCenter_(m_regionPSet.getParameter<double>("etaCenter")),
42  phiCenter_(m_regionPSet.getParameter<double>("phiCenter")),
43  deltaTrackEta(m_regionPSet.getParameter<double>("deltaEtaTrackRegion")),
44  deltaTrackPhi(m_regionPSet.getParameter<double>("deltaPhiTrackRegion")),
45  deltaL1JetEta(m_regionPSet.getParameter<double>("deltaEtaL1JetRegion")),
46  deltaL1JetPhi(m_regionPSet.getParameter<double>("deltaPhiL1JetRegion")),
47  usejets_(m_regionPSet.getParameter<bool>("useL1Jets")),
48  usetracks_(m_regionPSet.getParameter<bool>("useTracks")),
49  fixedReg_(m_regionPSet.getParameter<bool>("fixedReg")),
50  useIsoTracks_(m_regionPSet.getParameter<bool>("useIsoTracks")),
53  edm::LogVerbatim("HITRegionalPixelSeedGenerator") << "Enter the HITRegionalPixelSeedGenerator";
54 
55  if (usetracks_)
59  if (useIsoTracks_)
62  if (usejets_)
63  token_l1jet =
65  }
66 
67  ~HITRegionalPixelSeedGenerator() override = default;
68 
69  std::vector<std::unique_ptr<TrackingRegion> > regions(const edm::Event& e, const edm::EventSetup& es) const override {
70  std::vector<std::unique_ptr<TrackingRegion> > result;
71  float originz = 0.;
72 
73  double deltaZVertex = halflength;
74 
75  auto const& bfield = es.getData(token_bfield);
76  auto const& msmaker = es.getData(token_msmaker);
77 
78  if (usetracks_) {
80 
81  const reco::VertexCollection& vertCollection = e.get(token_vertex);
82  reco::VertexCollection::const_iterator ci = vertCollection.begin();
83 
84  if (!vertCollection.empty()) {
85  originz = ci->z();
86  } else {
87  deltaZVertex = 15.;
88  }
89 
90  GlobalVector globalVector(0, 0, 1);
91  if (tracks->empty())
92  return result;
93 
94  reco::TrackCollection::const_iterator itr = tracks->begin();
95  for (; itr != tracks->end(); itr++) {
96  GlobalVector ptrVec((itr)->px(), (itr)->py(), (itr)->pz());
97  globalVector = ptrVec;
98 
99  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(globalVector,
100  GlobalPoint(0, 0, originz),
101  ptmin,
102  originradius,
103  deltaZVertex,
106  bfield,
107  &msmaker));
108  }
109  }
110 
111  if (useIsoTracks_) {
113 
114  std::vector<edm::Ref<reco::IsolatedPixelTrackCandidateCollection> > isoPixTrackRefs;
115 
116  isotracks->getObjects(trigger::TriggerTrack, isoPixTrackRefs);
117 
118  const reco::VertexCollection& vertCollection = e.get(token_vertex);
119  reco::VertexCollection::const_iterator ci = vertCollection.begin();
120 
121  if (!vertCollection.empty()) {
122  originz = ci->z();
123  } else {
124  deltaZVertex = 15.;
125  }
126 
127  GlobalVector globalVector(0, 0, 1);
128  if (isoPixTrackRefs.empty())
129  return result;
130 
131  for (uint32_t p = 0; p < isoPixTrackRefs.size(); p++) {
132  GlobalVector ptrVec((isoPixTrackRefs[p]->track())->px(),
133  (isoPixTrackRefs[p]->track())->py(),
134  (isoPixTrackRefs[p]->track())->pz());
135  globalVector = ptrVec;
136 
137  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(globalVector,
138  GlobalPoint(0, 0, originz),
139  ptmin,
140  originradius,
141  deltaZVertex,
144  bfield,
145  &msmaker));
146  }
147  }
148 
149  if (usejets_) {
151  const reco::VertexCollection& vertCollection = e.get(token_vertex);
152  reco::VertexCollection::const_iterator ci = vertCollection.begin();
153  if (!vertCollection.empty()) {
154  originz = ci->z();
155  } else {
156  deltaZVertex = 15.;
157  }
158 
159  GlobalVector globalVector(0, 0, 1);
160  if (jets->empty())
161  return result;
162 
163  for (l1extra::L1JetParticleCollection::const_iterator iJet = jets->begin(); iJet != jets->end(); iJet++) {
164  GlobalVector jetVector(iJet->p4().x(), iJet->p4().y(), iJet->p4().z());
165  GlobalPoint vertex(0, 0, originz);
166 
167  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(
168  jetVector, vertex, ptmin, originradius, deltaZVertex, deltaL1JetEta, deltaL1JetPhi, bfield, &msmaker));
169  }
170  }
171  if (fixedReg_) {
172  GlobalVector fixedVector(cos(phiCenter_) * sin(2 * atan(exp(-etaCenter_))),
173  sin(phiCenter_) * sin(2 * atan(exp(-etaCenter_))),
174  cos(2 * atan(exp(-etaCenter_))));
175  GlobalPoint vertex(0, 0, originz);
176 
177  const reco::VertexCollection& vertCollection = e.get(token_vertex);
178  if (!vertCollection.empty()) {
179  // reco::VertexCollection::const_iterator ci = vertCollection.begin();
180  // originz = ci->z();
181  } else {
182  deltaZVertex = 15.;
183  }
184 
185  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(
186  fixedVector, vertex, ptmin, originradius, deltaZVertex, deltaL1JetEta, deltaL1JetPhi, bfield, &msmaker));
187  }
188 
189  return result;
190  }
191 
192 private:
194  const float ptmin;
195  const float originradius;
196  const float halflength;
197  const double etaCenter_;
198  const double phiCenter_;
199  const float deltaTrackEta;
200  const float deltaTrackPhi;
201  const float deltaL1JetEta;
202  const float deltaL1JetPhi;
203  const bool usejets_;
204  const bool usetracks_;
205  const bool fixedReg_;
206  const bool useIsoTracks_;
213 };
214 
219 #
Log< level::Info, true > LogVerbatim
edm::EDGetTokenT< reco::VertexCollection > token_vertex
edm::EDGetTokenT< l1extra::L1JetParticleCollection > token_l1jet
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > token_isoTrack
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
auto const & tracks
cannot be loose
Exp< T >::type exp(const T &t)
Definition: Exp.h:22
tuple result
Definition: mps_fire.py:311
bool getData(T &iHolder) const
Definition: EventSetup.h:122
Handle< PROD > getHandle(EDGetTokenT< PROD > token) const
Definition: Event.h:563
HITRegionalPixelSeedGenerator(const edm::ParameterSet &conf, edm::ConsumesCollector &&iC)
vector< PseudoJet > jets
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
std::vector< std::unique_ptr< TrackingRegion > > regions(const edm::Event &e, const edm::EventSetup &es) const override
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:346
~HITRegionalPixelSeedGenerator() override=default
const edm::ESGetToken< MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord > token_msmaker
std::vector< L1JetParticle > L1JetParticleCollection
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDGetTokenT< reco::TrackCollection > token_trks
#define DEFINE_EDM_PLUGIN(factory, type, name)
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > token_bfield
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283