CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HITRegionalPixelSeedGenerator.h
Go to the documentation of this file.
1 #ifndef HITRegionalPixelSeedGenerator_h
2 #define HITRegionalPixelSeedGenerator_h
3 
4 //
5 // Class: HITRegionalPixelSeedGenerator
6 
7 
25 
26 // Math
27 #include "Math/GenVector/VectorUtil.h"
28 #include "Math/GenVector/PxPyPzE4D.h"
29 
31 
34 
36  public:
37 
40  {
41  edm::LogInfo ("HITRegionalPixelSeedGenerator")<<"Enter the HITRegionalPixelSeedGenerator";
42 
43  edm::ParameterSet regionPSet = conf_.getParameter<edm::ParameterSet>("RegionPSet");
44 
45  ptmin=regionPSet.getParameter<double>("ptMin");
46  originradius=regionPSet.getParameter<double>("originRadius");
47  halflength=regionPSet.getParameter<double>("originHalfLength");
48  etaCenter_=regionPSet.getParameter<double>("etaCenter");
49  phiCenter_=regionPSet.getParameter<double>("phiCenter");
50  deltaTrackEta = regionPSet.getParameter<double>("deltaEtaTrackRegion");
51  deltaTrackPhi = regionPSet.getParameter<double>("deltaPhiTrackRegion");
52  deltaL1JetEta = regionPSet.getParameter<double>("deltaEtaL1JetRegion");
53  deltaL1JetPhi = regionPSet.getParameter<double>("deltaPhiL1JetRegion");
54  usejets_ = regionPSet.getParameter<bool>("useL1Jets");
55  usetracks_ = regionPSet.getParameter<bool>("useTracks");
56  useIsoTracks_ = regionPSet.getParameter<bool>("useIsoTracks");
57  fixedReg_ = regionPSet.getParameter<bool>("fixedReg");
58 
59  if (usetracks_) token_trks = iC.consumes<reco::TrackCollection>(regionPSet.getParameter<edm::InputTag>("trackSrc"));
60  if (usetracks_ || useIsoTracks_ || fixedReg_ || usejets_)
61  token_vertex = iC.consumes<reco::VertexCollection>(regionPSet.getParameter<edm::InputTag>("vertexSrc"));
62  if (useIsoTracks_) token_isoTrack = iC.consumes<trigger::TriggerFilterObjectWithRefs>(regionPSet.getParameter<edm::InputTag>("isoTrackSrc"));
63  if (usejets_) token_l1jet = iC.consumes<l1extra::L1JetParticleCollection>(regionPSet.getParameter<edm::InputTag>("l1tjetSrc"));
64  }
65 
67 
68 
69  virtual std::vector<std::unique_ptr<TrackingRegion> > regions(const edm::Event& e, const edm::EventSetup& es) const override
70  {
71  std::vector<std::unique_ptr<TrackingRegion> > result;
72  float originz =0.;
73 
74  double deltaZVertex = halflength;
75 
76 
77 
78  if (usetracks_)
79  {
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.size() > 0)
89  {
90  originz = ci->z();
91  }
92  else
93  {
94  deltaZVertex = 15.;
95  }
96 
97  GlobalVector globalVector(0,0,1);
98  if(tracks->size() == 0) return result;
99 
100  reco::TrackCollection::const_iterator itr = tracks->begin();
101  for(;itr != tracks->end();itr++)
102  {
103 
104  GlobalVector ptrVec((itr)->px(),(itr)->py(),(itr)->pz());
105  globalVector = ptrVec;
106 
107 
108  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(globalVector,
109  GlobalPoint(0,0,originz),
110  ptmin,
111  originradius,
112  deltaZVertex,
114  deltaTrackPhi));
115  }
116  }
117 
118  if (useIsoTracks_)
119  {
121  e.getByToken(token_isoTrack, isotracks);
122 
123  std::vector< edm::Ref<reco::IsolatedPixelTrackCandidateCollection> > isoPixTrackRefs;
124 
125  isotracks->getObjects(trigger::TriggerTrack, isoPixTrackRefs);
126 
128  e.getByToken(token_vertex,vertices);
129  const reco::VertexCollection vertCollection = *(vertices.product());
130  reco::VertexCollection::const_iterator ci = vertCollection.begin();
131 
132  if(vertCollection.size() > 0)
133  {
134  originz = ci->z();
135  }
136  else
137  {
138  deltaZVertex = 15.;
139  }
140 
141  GlobalVector globalVector(0,0,1);
142  if(isoPixTrackRefs.size() == 0) return result;
143 
144  for(uint32_t p=0; p<isoPixTrackRefs.size(); p++)
145  {
146  GlobalVector ptrVec((isoPixTrackRefs[p]->track())->px(),(isoPixTrackRefs[p]->track())->py(),(isoPixTrackRefs[p]->track())->pz());
147  globalVector = ptrVec;
148 
149  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(globalVector,
150  GlobalPoint(0,0,originz),
151  ptmin,
152  originradius,
153  deltaZVertex,
155  deltaTrackPhi));
156  }
157  }
158 
159  if (usejets_)
160  {
162  e.getByToken(token_l1jet, jets);
163 
165  e.getByToken(token_vertex, vertices);
166  const reco::VertexCollection vertCollection = *(vertices.product());
167  reco::VertexCollection::const_iterator ci = vertCollection.begin();
168  if(vertCollection.size() > 0)
169  {
170  originz = ci->z();
171  }
172  else
173  {
174  deltaZVertex = 15.;
175  }
176 
177  GlobalVector globalVector(0,0,1);
178  if(jets->size() == 0) return result;
179 
180  for (l1extra::L1JetParticleCollection::const_iterator iJet = jets->begin(); iJet != jets->end(); iJet++)
181  {
182  GlobalVector jetVector(iJet->p4().x(), iJet->p4().y(), iJet->p4().z());
183  GlobalPoint vertex(0, 0, originz);
184 
185  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>( jetVector,
186  vertex,
187  ptmin,
188  originradius,
189  deltaZVertex,
191  deltaL1JetPhi ));
192  }
193  }
194  if (fixedReg_)
195  {
196  GlobalVector fixedVector(cos(phiCenter_)*sin(2*atan(exp(-etaCenter_))), sin(phiCenter_)*sin(2*atan(exp(-etaCenter_))), cos(2*atan(exp(-etaCenter_))));
197  GlobalPoint vertex(0, 0, originz);
198 
200  e.getByToken(token_vertex,vertices);
201  const reco::VertexCollection vertCollection = *(vertices.product());
202  reco::VertexCollection::const_iterator ci = vertCollection.begin();
203  if(vertCollection.size() > 0)
204  {
205  originz = ci->z();
206  }
207  else
208  {
209  deltaZVertex = 15.;
210  }
211 
212  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>( fixedVector,
213  vertex,
214  ptmin,
215  originradius,
216  deltaZVertex,
218  deltaL1JetPhi ));
219  }
220 
221  return result;
222  }
223 
224 
225  private:
227 
228  float ptmin;
230  float halflength;
231  double etaCenter_;
232  double phiCenter_;
237  bool usejets_;
239  bool fixedReg_;
245 
246 };
247 
248 #endif
249 
250 
251 
T getParameter(std::string const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
edm::EDGetTokenT< reco::VertexCollection > token_vertex
edm::EDGetTokenT< l1extra::L1JetParticleCollection > token_l1jet
std::vector< L1JetParticle > L1JetParticleCollection
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
vector< PseudoJet > jets
tuple result
Definition: query.py:137
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
T const * product() const
Definition: Handle.h:81
tuple tracks
Definition: testEve_cfg.py:39
HITRegionalPixelSeedGenerator(const edm::ParameterSet &conf_, edm::ConsumesCollector &&iC)
edm::EDGetTokenT< reco::TrackCollection > token_trks
virtual std::vector< std::unique_ptr< TrackingRegion > > regions(const edm::Event &e, const edm::EventSetup &es) const override