CMS 3D CMS Logo

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