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 
25 
26 // Math
27 #include "Math/GenVector/VectorUtil.h"
28 #include "Math/GenVector/PxPyPzE4D.h"
29 
31 
34 
36 public:
39  edm::LogVerbatim("HITRegionalPixelSeedGenerator") << "Enter the HITRegionalPixelSeedGenerator";
40 
41  edm::ParameterSet regionPSet = conf.getParameter<edm::ParameterSet>("RegionPSet");
42 
43  ptmin = regionPSet.getParameter<double>("ptMin");
44  originradius = regionPSet.getParameter<double>("originRadius");
45  halflength = regionPSet.getParameter<double>("originHalfLength");
46  etaCenter_ = regionPSet.getParameter<double>("etaCenter");
47  phiCenter_ = regionPSet.getParameter<double>("phiCenter");
48  deltaTrackEta = regionPSet.getParameter<double>("deltaEtaTrackRegion");
49  deltaTrackPhi = regionPSet.getParameter<double>("deltaPhiTrackRegion");
50  deltaL1JetEta = regionPSet.getParameter<double>("deltaEtaL1JetRegion");
51  deltaL1JetPhi = regionPSet.getParameter<double>("deltaPhiL1JetRegion");
52  usejets_ = regionPSet.getParameter<bool>("useL1Jets");
53  usetracks_ = regionPSet.getParameter<bool>("useTracks");
54  useIsoTracks_ = regionPSet.getParameter<bool>("useIsoTracks");
55  fixedReg_ = regionPSet.getParameter<bool>("fixedReg");
56 
57  if (usetracks_)
58  token_trks = iC.consumes<reco::TrackCollection>(regionPSet.getParameter<edm::InputTag>("trackSrc"));
59  if (usetracks_ || useIsoTracks_ || fixedReg_ || usejets_)
60  token_vertex = iC.consumes<reco::VertexCollection>(regionPSet.getParameter<edm::InputTag>("vertexSrc"));
61  if (useIsoTracks_)
63  iC.consumes<trigger::TriggerFilterObjectWithRefs>(regionPSet.getParameter<edm::InputTag>("isoTrackSrc"));
64  if (usejets_)
65  token_l1jet = iC.consumes<l1extra::L1JetParticleCollection>(regionPSet.getParameter<edm::InputTag>("l1tjetSrc"));
66  }
67 
69 
70  std::vector<std::unique_ptr<TrackingRegion> > regions(const edm::Event& e, const edm::EventSetup& es) const override {
71  std::vector<std::unique_ptr<TrackingRegion> > result;
72  float originz = 0.;
73 
74  double deltaZVertex = halflength;
75 
76  auto const& bfield = es.getData(token_bfield);
77  auto const& msmaker = es.getData(token_msmaker);
78 
79  if (usetracks_) {
81  e.getByToken(token_trks, tracks);
82 
84  e.getByToken(token_vertex, vertices);
85  const reco::VertexCollection vertCollection = *(vertices.product());
86  reco::VertexCollection::const_iterator ci = vertCollection.begin();
87 
88  if (!vertCollection.empty()) {
89  originz = ci->z();
90  } else {
91  deltaZVertex = 15.;
92  }
93 
94  GlobalVector globalVector(0, 0, 1);
95  if (tracks->empty())
96  return result;
97 
98  reco::TrackCollection::const_iterator itr = tracks->begin();
99  for (; itr != tracks->end(); itr++) {
100  GlobalVector ptrVec((itr)->px(), (itr)->py(), (itr)->pz());
101  globalVector = ptrVec;
102 
103  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(globalVector,
104  GlobalPoint(0, 0, originz),
105  ptmin,
106  originradius,
107  deltaZVertex,
110  bfield,
111  &msmaker));
112  }
113  }
114 
115  if (useIsoTracks_) {
117  e.getByToken(token_isoTrack, isotracks);
118 
119  std::vector<edm::Ref<reco::IsolatedPixelTrackCandidateCollection> > isoPixTrackRefs;
120 
121  isotracks->getObjects(trigger::TriggerTrack, isoPixTrackRefs);
122 
124  e.getByToken(token_vertex, vertices);
125  const reco::VertexCollection vertCollection = *(vertices.product());
126  reco::VertexCollection::const_iterator ci = vertCollection.begin();
127 
128  if (!vertCollection.empty()) {
129  originz = ci->z();
130  } else {
131  deltaZVertex = 15.;
132  }
133 
134  GlobalVector globalVector(0, 0, 1);
135  if (isoPixTrackRefs.empty())
136  return result;
137 
138  for (uint32_t p = 0; p < isoPixTrackRefs.size(); p++) {
139  GlobalVector ptrVec((isoPixTrackRefs[p]->track())->px(),
140  (isoPixTrackRefs[p]->track())->py(),
141  (isoPixTrackRefs[p]->track())->pz());
142  globalVector = ptrVec;
143 
144  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(globalVector,
145  GlobalPoint(0, 0, originz),
146  ptmin,
147  originradius,
148  deltaZVertex,
151  bfield,
152  &msmaker));
153  }
154  }
155 
156  if (usejets_) {
158  e.getByToken(token_l1jet, jets);
159 
161  e.getByToken(token_vertex, vertices);
162  const reco::VertexCollection vertCollection = *(vertices.product());
163  reco::VertexCollection::const_iterator ci = vertCollection.begin();
164  if (!vertCollection.empty()) {
165  originz = ci->z();
166  } else {
167  deltaZVertex = 15.;
168  }
169 
170  GlobalVector globalVector(0, 0, 1);
171  if (jets->empty())
172  return result;
173 
174  for (l1extra::L1JetParticleCollection::const_iterator iJet = jets->begin(); iJet != jets->end(); iJet++) {
175  GlobalVector jetVector(iJet->p4().x(), iJet->p4().y(), iJet->p4().z());
176  GlobalPoint vertex(0, 0, originz);
177 
178  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(
179  jetVector, vertex, ptmin, originradius, deltaZVertex, deltaL1JetEta, deltaL1JetPhi, bfield, &msmaker));
180  }
181  }
182  if (fixedReg_) {
183  GlobalVector fixedVector(cos(phiCenter_) * sin(2 * atan(exp(-etaCenter_))),
184  sin(phiCenter_) * sin(2 * atan(exp(-etaCenter_))),
185  cos(2 * atan(exp(-etaCenter_))));
186  GlobalPoint vertex(0, 0, originz);
187 
189  e.getByToken(token_vertex, vertices);
190  const reco::VertexCollection vertCollection = *(vertices.product());
191  if (!vertCollection.empty()) {
192  // reco::VertexCollection::const_iterator ci = vertCollection.begin();
193  // originz = ci->z();
194  } else {
195  deltaZVertex = 15.;
196  }
197 
198  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(
199  fixedVector, vertex, ptmin, originradius, deltaZVertex, deltaL1JetEta, deltaL1JetPhi, bfield, &msmaker));
200  }
201 
202  return result;
203  }
204 
205 private:
206  float ptmin;
208  float halflength;
209  double etaCenter_;
210  double phiCenter_;
215  bool usejets_;
217  bool fixedReg_;
225 };
226 
231 #
Log< level::Info, true > LogVerbatim
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
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:128
edm::ESGetToken< MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord > token_msmaker
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
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > token_bfield
std::vector< L1JetParticle > L1JetParticleCollection
T const * product() const
Definition: Handle.h:70
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDGetTokenT< reco::TrackCollection > token_trks
#define DEFINE_EDM_PLUGIN(factory, type, name)
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283