CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
HITRegionalPixelSeedGenerator Class Reference

#include <HITRegionalPixelSeedGenerator.h>

Inheritance diagram for HITRegionalPixelSeedGenerator:
TrackingRegionProducer

Public Member Functions

 HITRegionalPixelSeedGenerator (const edm::ParameterSet &conf_, edm::ConsumesCollector &&iC)
 
virtual std::vector
< TrackingRegion * > 
regions (const edm::Event &e, const edm::EventSetup &es) const
 
virtual ~HITRegionalPixelSeedGenerator ()
 
- Public Member Functions inherited from TrackingRegionProducer
virtual ~TrackingRegionProducer ()
 

Private Attributes

edm::ParameterSet conf_
 
float deltaL1JetEta
 
float deltaL1JetPhi
 
float deltaTrackEta
 
float deltaTrackPhi
 
double etaCenter_
 
bool fixedReg_
 
float halflength
 
float originradius
 
double phiCenter_
 
float ptmin
 
edm::EDGetTokenT
< trigger::TriggerFilterObjectWithRefs
token_isoTrack
 
edm::EDGetTokenT
< l1extra::L1JetParticleCollection
token_l1jet
 
edm::EDGetTokenT
< reco::TrackCollection
token_trks
 
edm::EDGetTokenT
< reco::VertexCollection
token_vertex
 
bool useIsoTracks_
 
bool usejets_
 
bool usetracks_
 

Detailed Description

Definition at line 35 of file HITRegionalPixelSeedGenerator.h.

Constructor & Destructor Documentation

HITRegionalPixelSeedGenerator::HITRegionalPixelSeedGenerator ( const edm::ParameterSet conf_,
edm::ConsumesCollector &&  iC 
)
inlineexplicit

Definition at line 38 of file HITRegionalPixelSeedGenerator.h.

References deltaL1JetEta, deltaL1JetPhi, deltaTrackEta, deltaTrackPhi, etaCenter_, fixedReg_, edm::ParameterSet::getParameter(), halflength, originradius, phiCenter_, ptmin, token_isoTrack, token_l1jet, token_trks, token_vertex, useIsoTracks_, usejets_, and usetracks_.

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_)
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  }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
edm::EDGetTokenT< reco::VertexCollection > token_vertex
edm::EDGetTokenT< l1extra::L1JetParticleCollection > token_l1jet
std::vector< L1JetParticle > L1JetParticleCollection
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:13
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > token_isoTrack
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::EDGetTokenT< reco::TrackCollection > token_trks
virtual HITRegionalPixelSeedGenerator::~HITRegionalPixelSeedGenerator ( )
inlinevirtual

Definition at line 66 of file HITRegionalPixelSeedGenerator.h.

66 {}

Member Function Documentation

virtual std::vector<TrackingRegion* > HITRegionalPixelSeedGenerator::regions ( const edm::Event e,
const edm::EventSetup es 
) const
inlinevirtual

Implements TrackingRegionProducer.

Definition at line 69 of file HITRegionalPixelSeedGenerator.h.

References funct::cos(), deltaL1JetEta, deltaL1JetPhi, deltaTrackEta, deltaTrackPhi, etaCenter_, create_public_lumi_plots::exp, fixedReg_, edm::Event::getByToken(), halflength, fwrapper::jets, originradius, AlCaHLTBitMon_ParallelJobs::p, phiCenter_, edm::Handle< T >::product(), ptmin, query::result, funct::sin(), token_isoTrack, token_l1jet, token_trks, token_vertex, testEve_cfg::tracks, trigger::TriggerTrack, useIsoTracks_, usejets_, usetracks_, and HLT_25ns14e33_v1_cff::vertices.

70  {
71  std::vector<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  RectangularEtaPhiTrackingRegion* etaphiRegion = new RectangularEtaPhiTrackingRegion(globalVector,
109  GlobalPoint(0,0,originz),
110  ptmin,
111  originradius,
112  deltaZVertex,
114  deltaTrackPhi);
115  result.push_back(etaphiRegion);
116 
117  }
118  }
119 
120  if (useIsoTracks_)
121  {
123  e.getByToken(token_isoTrack, isotracks);
124 
125  std::vector< edm::Ref<reco::IsolatedPixelTrackCandidateCollection> > isoPixTrackRefs;
126 
127  isotracks->getObjects(trigger::TriggerTrack, isoPixTrackRefs);
128 
130  e.getByToken(token_vertex,vertices);
131  const reco::VertexCollection vertCollection = *(vertices.product());
132  reco::VertexCollection::const_iterator ci = vertCollection.begin();
133 
134  if(vertCollection.size() > 0)
135  {
136  originz = ci->z();
137  }
138  else
139  {
140  deltaZVertex = 15.;
141  }
142 
143  GlobalVector globalVector(0,0,1);
144  if(isoPixTrackRefs.size() == 0) return result;
145 
146  for(uint32_t p=0; p<isoPixTrackRefs.size(); p++)
147  {
148  GlobalVector ptrVec((isoPixTrackRefs[p]->track())->px(),(isoPixTrackRefs[p]->track())->py(),(isoPixTrackRefs[p]->track())->pz());
149  globalVector = ptrVec;
150 
151  RectangularEtaPhiTrackingRegion* etaphiRegion = new RectangularEtaPhiTrackingRegion(globalVector,
152  GlobalPoint(0,0,originz),
153  ptmin,
154  originradius,
155  deltaZVertex,
157  deltaTrackPhi);
158  result.push_back(etaphiRegion);
159  }
160  }
161 
162  if (usejets_)
163  {
165  e.getByToken(token_l1jet, jets);
166 
168  e.getByToken(token_vertex, vertices);
169  const reco::VertexCollection vertCollection = *(vertices.product());
170  reco::VertexCollection::const_iterator ci = vertCollection.begin();
171  if(vertCollection.size() > 0)
172  {
173  originz = ci->z();
174  }
175  else
176  {
177  deltaZVertex = 15.;
178  }
179 
180  GlobalVector globalVector(0,0,1);
181  if(jets->size() == 0) return result;
182 
183  for (l1extra::L1JetParticleCollection::const_iterator iJet = jets->begin(); iJet != jets->end(); iJet++)
184  {
185  GlobalVector jetVector(iJet->p4().x(), iJet->p4().y(), iJet->p4().z());
186  GlobalPoint vertex(0, 0, originz);
187 
189  vertex,
190  ptmin,
191  originradius,
192  deltaZVertex,
194  deltaL1JetPhi );
195  result.push_back(etaphiRegion);
196  }
197  }
198  if (fixedReg_)
199  {
200  GlobalVector fixedVector(cos(phiCenter_)*sin(2*atan(exp(-etaCenter_))), sin(phiCenter_)*sin(2*atan(exp(-etaCenter_))), cos(2*atan(exp(-etaCenter_))));
201  GlobalPoint vertex(0, 0, originz);
202 
204  e.getByToken(token_vertex,vertices);
205  const reco::VertexCollection vertCollection = *(vertices.product());
206  reco::VertexCollection::const_iterator ci = vertCollection.begin();
207  if(vertCollection.size() > 0)
208  {
209  originz = ci->z();
210  }
211  else
212  {
213  deltaZVertex = 15.;
214  }
215 
216  RectangularEtaPhiTrackingRegion* etaphiRegion = new RectangularEtaPhiTrackingRegion( fixedVector,
217  vertex,
218  ptmin,
219  originradius,
220  deltaZVertex,
222  deltaL1JetPhi );
223  result.push_back(etaphiRegion);
224  }
225 
226  return result;
227  }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
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
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
edm::EDGetTokenT< reco::TrackCollection > token_trks

Member Data Documentation

edm::ParameterSet HITRegionalPixelSeedGenerator::conf_
private

Definition at line 231 of file HITRegionalPixelSeedGenerator.h.

float HITRegionalPixelSeedGenerator::deltaL1JetEta
private

Definition at line 240 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

float HITRegionalPixelSeedGenerator::deltaL1JetPhi
private

Definition at line 241 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

float HITRegionalPixelSeedGenerator::deltaTrackEta
private

Definition at line 238 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

float HITRegionalPixelSeedGenerator::deltaTrackPhi
private

Definition at line 239 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

double HITRegionalPixelSeedGenerator::etaCenter_
private

Definition at line 236 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

bool HITRegionalPixelSeedGenerator::fixedReg_
private

Definition at line 244 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

float HITRegionalPixelSeedGenerator::halflength
private

Definition at line 235 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

float HITRegionalPixelSeedGenerator::originradius
private

Definition at line 234 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

double HITRegionalPixelSeedGenerator::phiCenter_
private

Definition at line 237 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

float HITRegionalPixelSeedGenerator::ptmin
private

Definition at line 233 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

edm::EDGetTokenT<trigger::TriggerFilterObjectWithRefs> HITRegionalPixelSeedGenerator::token_isoTrack
private

Definition at line 248 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

edm::EDGetTokenT<l1extra::L1JetParticleCollection> HITRegionalPixelSeedGenerator::token_l1jet
private

Definition at line 249 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

edm::EDGetTokenT<reco::TrackCollection> HITRegionalPixelSeedGenerator::token_trks
private

Definition at line 246 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

edm::EDGetTokenT<reco::VertexCollection> HITRegionalPixelSeedGenerator::token_vertex
private

Definition at line 247 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

bool HITRegionalPixelSeedGenerator::useIsoTracks_
private

Definition at line 245 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

bool HITRegionalPixelSeedGenerator::usejets_
private

Definition at line 242 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

bool HITRegionalPixelSeedGenerator::usetracks_
private

Definition at line 243 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().