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_)
 
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
 
edm::InputTag isoTracksrc_
 
edm::InputTag l1jetsrc_
 
float originradius
 
double phiCenter_
 
float ptmin
 
edm::InputTag tracksrc_
 
bool useIsoTracks_
 
bool usejets_
 
bool usetracks_
 
std::string vertexSrc
 

Detailed Description

Definition at line 34 of file HITRegionalPixelSeedGenerator.h.

Constructor & Destructor Documentation

HITRegionalPixelSeedGenerator::HITRegionalPixelSeedGenerator ( const edm::ParameterSet conf_)
inlineexplicit

Definition at line 37 of file HITRegionalPixelSeedGenerator.h.

References deltaL1JetEta, deltaL1JetPhi, deltaTrackEta, deltaTrackPhi, etaCenter_, fixedReg_, edm::ParameterSet::getParameter(), halflength, isoTracksrc_, l1jetsrc_, originradius, phiCenter_, ptmin, tracksrc_, useIsoTracks_, usejets_, usetracks_, and vertexSrc.

38  {
39  edm::LogInfo ("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  vertexSrc=regionPSet.getParameter<std::string>("vertexSrc");
47  etaCenter_=regionPSet.getParameter<double>("etaCenter");
48  phiCenter_=regionPSet.getParameter<double>("phiCenter");
49  deltaTrackEta = regionPSet.getParameter<double>("deltaEtaTrackRegion");
50  deltaTrackPhi = regionPSet.getParameter<double>("deltaPhiTrackRegion");
51  deltaL1JetEta = regionPSet.getParameter<double>("deltaEtaL1JetRegion");
52  deltaL1JetPhi = regionPSet.getParameter<double>("deltaPhiL1JetRegion");
53  tracksrc_ = regionPSet.getParameter<edm::InputTag>("trackSrc");
54  isoTracksrc_ = regionPSet.getParameter<edm::InputTag>("isoTrackSrc");
55  l1jetsrc_ = regionPSet.getParameter<edm::InputTag>("l1tjetSrc");
56  usejets_ = regionPSet.getParameter<bool>("useL1Jets");
57  usetracks_ = regionPSet.getParameter<bool>("useTracks");
58  useIsoTracks_ = regionPSet.getParameter<bool>("useIsoTracks");
59  fixedReg_ = regionPSet.getParameter<bool>("fixedReg");
60  }
T getParameter(std::string const &) const
virtual HITRegionalPixelSeedGenerator::~HITRegionalPixelSeedGenerator ( )
inlinevirtual

Definition at line 62 of file HITRegionalPixelSeedGenerator.h.

62 {}

Member Function Documentation

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

Implements TrackingRegionProducer.

Definition at line 65 of file HITRegionalPixelSeedGenerator.h.

References funct::cos(), deltaL1JetEta, deltaL1JetPhi, deltaTrackEta, deltaTrackPhi, etaCenter_, create_public_lumi_plots::exp, fixedReg_, edm::Event::getByLabel(), halflength, isoTracksrc_, fwrapper::jets, l1jetsrc_, originradius, AlCaHLTBitMon_ParallelJobs::p, phiCenter_, edm::Handle< T >::product(), ptmin, query::result, funct::sin(), testEve_cfg::tracks, tracksrc_, trigger::TriggerTrack, useIsoTracks_, usejets_, usetracks_, and vertexSrc.

66  {
67  std::vector<TrackingRegion* > result;
68  float originz =0.;
69 
70  double deltaZVertex = halflength;
71 
72 
73 
74  if (usetracks_)
75  {
77  e.getByLabel(tracksrc_, tracks);
78 
80  e.getByLabel(vertexSrc,vertices);
81  const reco::VertexCollection vertCollection = *(vertices.product());
82  reco::VertexCollection::const_iterator ci = vertCollection.begin();
83 
84  if(vertCollection.size() > 0)
85  {
86  originz = ci->z();
87  }
88  else
89  {
90  deltaZVertex = 15.;
91  }
92 
93  GlobalVector globalVector(0,0,1);
94  if(tracks->size() == 0) return result;
95 
96  reco::TrackCollection::const_iterator itr = tracks->begin();
97  for(;itr != tracks->end();itr++)
98  {
99 
100  GlobalVector ptrVec((itr)->px(),(itr)->py(),(itr)->pz());
101  globalVector = ptrVec;
102 
103 
104  RectangularEtaPhiTrackingRegion* etaphiRegion = new RectangularEtaPhiTrackingRegion(globalVector,
105  GlobalPoint(0,0,originz),
106  ptmin,
107  originradius,
108  deltaZVertex,
110  deltaTrackPhi);
111  result.push_back(etaphiRegion);
112 
113  }
114  }
115 
116  if (useIsoTracks_)
117  {
119  e.getByLabel(isoTracksrc_, isotracks);
120 
121  std::vector< edm::Ref<reco::IsolatedPixelTrackCandidateCollection> > isoPixTrackRefs;
122 
123  isotracks->getObjects(trigger::TriggerTrack, isoPixTrackRefs);
124 
126  e.getByLabel(vertexSrc,vertices);
127  const reco::VertexCollection vertCollection = *(vertices.product());
128  reco::VertexCollection::const_iterator ci = vertCollection.begin();
129 
130  if(vertCollection.size() > 0)
131  {
132  originz = ci->z();
133  }
134  else
135  {
136  deltaZVertex = 15.;
137  }
138 
139  GlobalVector globalVector(0,0,1);
140  if(isoPixTrackRefs.size() == 0) return result;
141 
142  for(uint32_t p=0; p<isoPixTrackRefs.size(); p++)
143  {
144  GlobalVector ptrVec((isoPixTrackRefs[p]->track())->px(),(isoPixTrackRefs[p]->track())->py(),(isoPixTrackRefs[p]->track())->pz());
145  globalVector = ptrVec;
146 
147  RectangularEtaPhiTrackingRegion* etaphiRegion = new RectangularEtaPhiTrackingRegion(globalVector,
148  GlobalPoint(0,0,originz),
149  ptmin,
150  originradius,
151  deltaZVertex,
153  deltaTrackPhi);
154  result.push_back(etaphiRegion);
155  }
156  }
157 
158  if (usejets_)
159  {
161  e.getByLabel(l1jetsrc_, jets);
162 
164  e.getByLabel(vertexSrc,vertices);
165  const reco::VertexCollection vertCollection = *(vertices.product());
166  reco::VertexCollection::const_iterator ci = vertCollection.begin();
167  if(vertCollection.size() > 0)
168  {
169  originz = ci->z();
170  }
171  else
172  {
173  deltaZVertex = 15.;
174  }
175 
176  GlobalVector globalVector(0,0,1);
177  if(jets->size() == 0) return result;
178 
179  for (l1extra::L1JetParticleCollection::const_iterator iJet = jets->begin(); iJet != jets->end(); iJet++)
180  {
181  GlobalVector jetVector(iJet->p4().x(), iJet->p4().y(), iJet->p4().z());
182  GlobalPoint vertex(0, 0, originz);
183 
185  vertex,
186  ptmin,
187  originradius,
188  deltaZVertex,
190  deltaL1JetPhi );
191  result.push_back(etaphiRegion);
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.getByLabel(vertexSrc,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  RectangularEtaPhiTrackingRegion* etaphiRegion = new RectangularEtaPhiTrackingRegion( fixedVector,
213  vertex,
214  ptmin,
215  originradius,
216  deltaZVertex,
218  deltaL1JetPhi );
219  result.push_back(etaphiRegion);
220  }
221 
222  return result;
223  }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
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
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
tuple tracks
Definition: testEve_cfg.py:39
T const * product() const
Definition: Handle.h:74

Member Data Documentation

edm::ParameterSet HITRegionalPixelSeedGenerator::conf_
private

Definition at line 227 of file HITRegionalPixelSeedGenerator.h.

float HITRegionalPixelSeedGenerator::deltaL1JetEta
private

Definition at line 236 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

float HITRegionalPixelSeedGenerator::deltaL1JetPhi
private

Definition at line 237 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

float HITRegionalPixelSeedGenerator::deltaTrackEta
private

Definition at line 234 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

float HITRegionalPixelSeedGenerator::deltaTrackPhi
private

Definition at line 235 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

double HITRegionalPixelSeedGenerator::etaCenter_
private

Definition at line 232 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 231 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

edm::InputTag HITRegionalPixelSeedGenerator::isoTracksrc_
private

Definition at line 239 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

edm::InputTag HITRegionalPixelSeedGenerator::l1jetsrc_
private

Definition at line 241 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

float HITRegionalPixelSeedGenerator::originradius
private

Definition at line 230 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

double HITRegionalPixelSeedGenerator::phiCenter_
private

Definition at line 233 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

float HITRegionalPixelSeedGenerator::ptmin
private

Definition at line 229 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

edm::InputTag HITRegionalPixelSeedGenerator::tracksrc_
private

Definition at line 238 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().

std::string HITRegionalPixelSeedGenerator::vertexSrc
private

Definition at line 240 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().