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 35 of file HITRegionalPixelSeedGenerator.h.

Constructor & Destructor Documentation

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

Definition at line 38 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.

39  {
40  edm::LogInfo ("HITRegionalPixelSeedGenerator")<<"Enter the HITRegionalPixelSeedGenerator";
41 
42  edm::ParameterSet regionPSet = conf_.getParameter<edm::ParameterSet>("RegionPSet");
43 
44  ptmin=regionPSet.getParameter<double>("ptMin");
45  originradius=regionPSet.getParameter<double>("originRadius");
46  halflength=regionPSet.getParameter<double>("originHalfLength");
47  vertexSrc=regionPSet.getParameter<std::string>("vertexSrc");
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  tracksrc_ = regionPSet.getParameter<edm::InputTag>("trackSrc");
55  isoTracksrc_ = regionPSet.getParameter<edm::InputTag>("isoTrackSrc");
56  l1jetsrc_ = regionPSet.getParameter<edm::InputTag>("l1tjetSrc");
57  usejets_ = regionPSet.getParameter<bool>("useL1Jets");
58  usetracks_ = regionPSet.getParameter<bool>("useTracks");
59  useIsoTracks_ = regionPSet.getParameter<bool>("useIsoTracks");
60  fixedReg_ = regionPSet.getParameter<bool>("fixedReg");
61  }
T getParameter(std::string const &) const
virtual HITRegionalPixelSeedGenerator::~HITRegionalPixelSeedGenerator ( )
inlinevirtual

Definition at line 63 of file HITRegionalPixelSeedGenerator.h.

63 {}

Member Function Documentation

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

Implements TrackingRegionProducer.

Definition at line 66 of file HITRegionalPixelSeedGenerator.h.

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

67  {
68  std::vector<TrackingRegion* > result;
69  float originz =0.;
70 
71  double deltaZVertex = halflength;
72 
73 
74 
75  if (usetracks_)
76  {
78  e.getByLabel(tracksrc_, tracks);
79 
81  e.getByLabel(vertexSrc,vertices);
82  const reco::VertexCollection vertCollection = *(vertices.product());
83  reco::VertexCollection::const_iterator ci = vertCollection.begin();
84 
85  if(vertCollection.size() > 0)
86  {
87  originz = ci->z();
88  }
89  else
90  {
91  deltaZVertex = 15.;
92  }
93 
94  GlobalVector globalVector(0,0,1);
95  if(tracks->size() == 0) return result;
96 
97  reco::TrackCollection::const_iterator itr = tracks->begin();
98  for(;itr != tracks->end();itr++)
99  {
100 
101  GlobalVector ptrVec((itr)->px(),(itr)->py(),(itr)->pz());
102  globalVector = ptrVec;
103 
104 
105  RectangularEtaPhiTrackingRegion* etaphiRegion = new RectangularEtaPhiTrackingRegion(globalVector,
106  GlobalPoint(0,0,originz),
107  ptmin,
108  originradius,
109  deltaZVertex,
111  deltaTrackPhi);
112  result.push_back(etaphiRegion);
113 
114  }
115  }
116 
117  if (useIsoTracks_)
118  {
120  e.getByLabel(isoTracksrc_, isotracks);
121 
122  std::vector< edm::Ref<reco::IsolatedPixelTrackCandidateCollection> > isoPixTrackRefs;
123 
124  isotracks->getObjects(trigger::TriggerTrack, isoPixTrackRefs);
125 
127  e.getByLabel(vertexSrc,vertices);
128  const reco::VertexCollection vertCollection = *(vertices.product());
129  reco::VertexCollection::const_iterator ci = vertCollection.begin();
130 
131  if(vertCollection.size() > 0)
132  {
133  originz = ci->z();
134  }
135  else
136  {
137  deltaZVertex = 15.;
138  }
139 
140  GlobalVector globalVector(0,0,1);
141  if(isoPixTrackRefs.size() == 0) return result;
142 
143  for(uint32_t p=0; p<isoPixTrackRefs.size(); p++)
144  {
145  GlobalVector ptrVec((isoPixTrackRefs[p]->track())->px(),(isoPixTrackRefs[p]->track())->py(),(isoPixTrackRefs[p]->track())->pz());
146  globalVector = ptrVec;
147 
148  RectangularEtaPhiTrackingRegion* etaphiRegion = new RectangularEtaPhiTrackingRegion(globalVector,
149  GlobalPoint(0,0,originz),
150  ptmin,
151  originradius,
152  deltaZVertex,
154  deltaTrackPhi);
155  result.push_back(etaphiRegion);
156  }
157  }
158 
159  if (usejets_)
160  {
162  e.getByLabel(l1jetsrc_, jets);
163 
165  e.getByLabel(vertexSrc,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 
186  vertex,
187  ptmin,
188  originradius,
189  deltaZVertex,
191  deltaL1JetPhi );
192  result.push_back(etaphiRegion);
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.getByLabel(vertexSrc,vertices);
202  const reco::VertexCollection vertCollection = *(vertices.product());
203  reco::VertexCollection::const_iterator ci = vertCollection.begin();
204  if(vertCollection.size() > 0)
205  {
206  originz = ci->z();
207  }
208  else
209  {
210  deltaZVertex = 15.;
211  }
212 
213  RectangularEtaPhiTrackingRegion* etaphiRegion = new RectangularEtaPhiTrackingRegion( fixedVector,
214  vertex,
215  ptmin,
216  originradius,
217  deltaZVertex,
219  deltaL1JetPhi );
220  result.push_back(etaphiRegion);
221  }
222 
223  return result;
224  }
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
Exp< T >::type exp(const T &t)
Definition: Exp.h:22
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:359
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 228 of file HITRegionalPixelSeedGenerator.h.

float HITRegionalPixelSeedGenerator::deltaL1JetEta
private

Definition at line 237 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

float HITRegionalPixelSeedGenerator::deltaL1JetPhi
private

Definition at line 238 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

float HITRegionalPixelSeedGenerator::deltaTrackEta
private

Definition at line 235 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

float HITRegionalPixelSeedGenerator::deltaTrackPhi
private

Definition at line 236 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

double HITRegionalPixelSeedGenerator::etaCenter_
private

Definition at line 233 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

bool HITRegionalPixelSeedGenerator::fixedReg_
private

Definition at line 245 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

float HITRegionalPixelSeedGenerator::halflength
private

Definition at line 232 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

edm::InputTag HITRegionalPixelSeedGenerator::isoTracksrc_
private

Definition at line 240 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

edm::InputTag HITRegionalPixelSeedGenerator::l1jetsrc_
private

Definition at line 242 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

float HITRegionalPixelSeedGenerator::originradius
private

Definition at line 231 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

double HITRegionalPixelSeedGenerator::phiCenter_
private

Definition at line 234 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

float HITRegionalPixelSeedGenerator::ptmin
private

Definition at line 230 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

edm::InputTag HITRegionalPixelSeedGenerator::tracksrc_
private

Definition at line 239 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

bool HITRegionalPixelSeedGenerator::useIsoTracks_
private

Definition at line 246 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

bool HITRegionalPixelSeedGenerator::usejets_
private

Definition at line 243 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

bool HITRegionalPixelSeedGenerator::usetracks_
private

Definition at line 244 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().

std::string HITRegionalPixelSeedGenerator::vertexSrc
private

Definition at line 241 of file HITRegionalPixelSeedGenerator.h.

Referenced by HITRegionalPixelSeedGenerator(), and regions().