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
TauRegionalPixelSeedGenerator Class Reference

#include <TauRegionalPixelSeedGenerator.h>

Inheritance diagram for TauRegionalPixelSeedGenerator:
TrackingRegionProducer

Public Member Functions

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

Private Attributes

edm::ParameterSet conf_
 
float m_deltaEta
 
float m_deltaPhi
 
float m_halfLength
 
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_howToUseMeasurementTracker
 
float m_originRadius
 
float m_ptMin
 
bool m_searchOpt
 
edm::EDGetTokenT
< reco::CandidateView
token_jet
 
edm::EDGetTokenT
< MeasurementTrackerEvent
token_measurementTracker
 
edm::EDGetTokenT
< reco::VertexCollection
token_vertex
 

Detailed Description

Definition at line 33 of file TauRegionalPixelSeedGenerator.h.

Constructor & Destructor Documentation

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

Definition at line 36 of file TauRegionalPixelSeedGenerator.h.

References edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), RectangularEtaPhiTrackingRegion::kNever, m_deltaEta, m_deltaPhi, m_halfLength, m_howToUseMeasurementTracker, m_originRadius, m_ptMin, m_searchOpt, AlCaHLTBitMon_QueryRunRegistry::string, RectangularEtaPhiTrackingRegion::stringToUseMeasurementTracker(), token_jet, token_measurementTracker, and token_vertex.

37  {
38  edm::LogInfo ("TauRegionalPixelSeedGenerator") << "Enter the TauRegionalPixelSeedGenerator";
39 
40  edm::ParameterSet regionPSet = conf_.getParameter<edm::ParameterSet>("RegionPSet");
41 
42  m_ptMin = regionPSet.getParameter<double>("ptMin");
43  m_originRadius = regionPSet.getParameter<double>("originRadius");
44  m_halfLength = regionPSet.getParameter<double>("originHalfLength");
45  m_deltaEta = regionPSet.getParameter<double>("deltaEtaRegion");
46  m_deltaPhi = regionPSet.getParameter<double>("deltaPhiRegion");
48  token_vertex = iC.consumes<reco::VertexCollection>(regionPSet.getParameter<edm::InputTag>("vertexSrc"));
49  if (regionPSet.exists("searchOpt")){
50  m_searchOpt = regionPSet.getParameter<bool>("searchOpt");
51  }
52  else{
53  m_searchOpt = false;
54  }
57  token_measurementTracker = iC.consumes<MeasurementTrackerEvent>(regionPSet.getParameter<std::string>("measurementTrackerName"));
58  }
59  }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
static UseMeasurementTracker stringToUseMeasurementTracker(const std::string &name)
edm::EDGetTokenT< MeasurementTrackerEvent > token_measurementTracker
edm::EDGetTokenT< reco::VertexCollection > token_vertex
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_howToUseMeasurementTracker
edm::EDGetTokenT< reco::CandidateView > token_jet
virtual TauRegionalPixelSeedGenerator::~TauRegionalPixelSeedGenerator ( )
inlinevirtual

Definition at line 61 of file TauRegionalPixelSeedGenerator.h.

61 {}

Member Function Documentation

virtual std::vector<std::unique_ptr<TrackingRegion> > TauRegionalPixelSeedGenerator::regions ( const edm::Event e,
const edm::EventSetup es 
) const
inlineoverridevirtual

Implements TrackingRegionProducer.

Definition at line 64 of file TauRegionalPixelSeedGenerator.h.

References edm::Event::getByToken(), edm::EDGetTokenT< T >::isUninitialized(), m_deltaEta, m_deltaPhi, m_halfLength, m_howToUseMeasurementTracker, m_originRadius, m_ptMin, m_searchOpt, HLT_FULL_cff::measurementTracker, edm::Handle< T >::product(), query::result, token_jet, token_measurementTracker, token_vertex, and HLT_FULL_cff::vertices.

64  {
65  std::vector<std::unique_ptr<TrackingRegion> > result;
66 
67  // double originZ;
68  double deltaZVertex, deltaRho;
69  GlobalPoint vertex;
70  // get the primary vertex
72  e.getByToken(token_vertex, h_vertices);
73  const reco::VertexCollection & vertices = * h_vertices;
74  if (not vertices.empty()) {
75 // originZ = vertices.front().z();
76  GlobalPoint myTmp(vertices.at(0).position().x(),vertices.at(0).position().y(), vertices.at(0).position().z());
77  vertex = myTmp;
78  deltaZVertex = m_halfLength;
79  deltaRho = m_originRadius;
80  } else {
81  // originZ = 0.;
82  GlobalPoint myTmp(0.,0.,0.);
83  vertex = myTmp;
84  deltaRho = 1.;
85  deltaZVertex = 15.;
86  }
87 
88  // get the jet direction
90  e.getByToken(token_jet, h_jets);
91 
96  measurementTracker = hmte.product();
97  }
98 
99  for(const reco::Candidate& myJet: *h_jets)
100  {
101  GlobalVector jetVector(myJet.momentum().x(),myJet.momentum().y(),myJet.momentum().z());
102 // GlobalPoint vertex(0, 0, originZ);
103  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>( jetVector,
104  vertex,
105  m_ptMin,
106  deltaRho,
107  deltaZVertex,
108  m_deltaEta,
109  m_deltaPhi,
111  true,
112  measurementTracker,
113  m_searchOpt));
114  }
115 
116  return result;
117  }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
tuple measurementTracker
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::EDGetTokenT< MeasurementTrackerEvent > token_measurementTracker
tuple result
Definition: query.py:137
edm::EDGetTokenT< reco::VertexCollection > token_vertex
T const * product() const
Definition: Handle.h:81
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_howToUseMeasurementTracker
bool isUninitialized() const
Definition: EDGetToken.h:73
edm::EDGetTokenT< reco::CandidateView > token_jet

Member Data Documentation

edm::ParameterSet TauRegionalPixelSeedGenerator::conf_
private

Definition at line 120 of file TauRegionalPixelSeedGenerator.h.

float TauRegionalPixelSeedGenerator::m_deltaEta
private

Definition at line 125 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

float TauRegionalPixelSeedGenerator::m_deltaPhi
private

Definition at line 126 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

float TauRegionalPixelSeedGenerator::m_halfLength
private

Definition at line 124 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

RectangularEtaPhiTrackingRegion::UseMeasurementTracker TauRegionalPixelSeedGenerator::m_howToUseMeasurementTracker
private

Definition at line 130 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

float TauRegionalPixelSeedGenerator::m_originRadius
private

Definition at line 123 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

float TauRegionalPixelSeedGenerator::m_ptMin
private

Definition at line 122 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

bool TauRegionalPixelSeedGenerator::m_searchOpt
private

Definition at line 131 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

edm::EDGetTokenT<reco::CandidateView> TauRegionalPixelSeedGenerator::token_jet
private

Definition at line 128 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

edm::EDGetTokenT<MeasurementTrackerEvent> TauRegionalPixelSeedGenerator::token_measurementTracker
private

Definition at line 129 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

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

Definition at line 127 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().