CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Attributes
TauRegionalPixelSeedGenerator Class Reference

#include <TauRegionalPixelSeedGenerator.h>

Inheritance diagram for TauRegionalPixelSeedGenerator:
TrackingRegionProducer

Public Member Functions

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)
 
 ~TauRegionalPixelSeedGenerator () override
 
- Public Member Functions inherited from TrackingRegionProducer
virtual ~TrackingRegionProducer ()
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Attributes

float m_deltaEta
 
float m_deltaPhi
 
float m_halfLength
 
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_howToUseMeasurementTracker
 
float m_originRadius
 
float m_ptMin
 
bool m_searchOpt
 
edm::ESGetToken< MagneticField, IdealMagneticFieldRecordtoken_bfield
 
edm::EDGetTokenT< reco::CandidateViewtoken_jet
 
edm::EDGetTokenT< MeasurementTrackerEventtoken_measurementTracker
 
edm::ESGetToken< MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecordtoken_msmaker
 
edm::EDGetTokenT< reco::VertexCollectiontoken_vertex
 

Detailed Description

Definition at line 36 of file TauRegionalPixelSeedGenerator.h.

Constructor & Destructor Documentation

◆ TauRegionalPixelSeedGenerator()

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

Definition at line 38 of file TauRegionalPixelSeedGenerator.h.

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

40  edm::LogInfo("TauRegionalPixelSeedGenerator") << "Enter the TauRegionalPixelSeedGenerator";
41 
42  edm::ParameterSet regionPSet = conf.getParameter<edm::ParameterSet>("RegionPSet");
43 
44  m_ptMin = regionPSet.getParameter<double>("ptMin");
45  m_originRadius = regionPSet.getParameter<double>("originRadius");
46  m_halfLength = regionPSet.getParameter<double>("originHalfLength");
47  m_deltaEta = regionPSet.getParameter<double>("deltaEtaRegion");
48  m_deltaPhi = regionPSet.getParameter<double>("deltaPhiRegion");
50  token_vertex = iC.consumes<reco::VertexCollection>(regionPSet.getParameter<edm::InputTag>("vertexSrc"));
51  if (regionPSet.exists("searchOpt")) {
52  m_searchOpt = regionPSet.getParameter<bool>("searchOpt");
53  } else {
54  m_searchOpt = false;
55  }
57  regionPSet.getParameter<std::string>("howToUseMeasurementTracker"));
59  // temporary until everything migrated to InputTag
61  if (regionPSet.existsAs<edm::InputTag>("measurementTrackerName")) {
62  tag = regionPSet.getParameter<edm::InputTag>("measurementTrackerName");
63  } else {
64  tag = edm::InputTag(regionPSet.getParameter<std::string>("measurementTrackerName"));
65  }
66 
68  }
69  }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
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::ESGetToken< MagneticField, IdealMagneticFieldRecord > token_bfield
edm::EDGetTokenT< reco::VertexCollection > token_vertex
Log< level::Info, false > LogInfo
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_howToUseMeasurementTracker
edm::ESGetToken< MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord > token_msmaker
edm::EDGetTokenT< reco::CandidateView > token_jet

◆ ~TauRegionalPixelSeedGenerator()

TauRegionalPixelSeedGenerator::~TauRegionalPixelSeedGenerator ( )
inlineoverride

Definition at line 71 of file TauRegionalPixelSeedGenerator.h.

71 {}

Member Function Documentation

◆ fillDescriptions()

static void TauRegionalPixelSeedGenerator::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
inlinestatic

Definition at line 73 of file TauRegionalPixelSeedGenerator.h.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), submitPVResolutionJobs::desc, ProducerED_cfi::InputTag, and AlCaHLTBitMon_QueryRunRegistry::string.

73  {
75 
76  desc.add<double>("ptMin", 5.0);
77  desc.add<double>("originHalfLength", 0.2);
78  desc.add<double>("originRadius", 0.2);
79  desc.add<double>("deltaEtaRegion", 0.1);
80  desc.add<double>("deltaPhiRegion", 0.1);
81  desc.add<edm::InputTag>("JetSrc", edm::InputTag("icone5Tau1"));
82  desc.add<edm::InputTag>("vertexSrc", edm::InputTag("pixelVertices"));
83  desc.add<bool>("searchOpt", false);
84 
85  desc.add<std::string>("howToUseMeasurementTracker", "ForSiStrips");
86 
87  // allow both InputTag and string for the moment, use InputTag as the default
89  "measurementTrackerName", edm::InputTag("MeasurementTrackerEvent"), true) xor
90  edm::ParameterDescription<std::string>("measurementTrackerName", "MeasurementTrackerEvent", true));
91 
92  // Only for backwards-compatibility
94  descRegion.add<edm::ParameterSetDescription>("RegionPSet", desc);
95 
96  descriptions.add("tauRegionalPixelSeedTrackingRegions", descRegion);
97  }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)

◆ regions()

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

Implements TrackingRegionProducer.

Definition at line 99 of file TauRegionalPixelSeedGenerator.h.

References MillePedeFileConverter_cfg::e, edm::EventSetup::getData(), edm::EDGetTokenT< T >::isUninitialized(), m_deltaEta, m_deltaPhi, m_halfLength, m_howToUseMeasurementTracker, m_originRadius, m_ptMin, m_searchOpt, HLTSiStripMonitoring_cff::measurementTracker, edm::Handle< T >::product(), mps_fire::result, token_bfield, token_jet, token_measurementTracker, token_msmaker, token_vertex, bphysicsOniaDQM_cfi::vertex, and AlignmentTracksFromVertexSelector_cfi::vertices.

99  {
100  std::vector<std::unique_ptr<TrackingRegion> > result;
101 
102  // double originZ;
103  double deltaZVertex, deltaRho;
105  // get the primary vertex
107  e.getByToken(token_vertex, h_vertices);
108  const reco::VertexCollection& vertices = *h_vertices;
109  if (not vertices.empty()) {
110  // originZ = vertices.front().z();
111  GlobalPoint myTmp(vertices.at(0).position().x(), vertices.at(0).position().y(), vertices.at(0).position().z());
112  vertex = myTmp;
113  deltaZVertex = m_halfLength;
114  deltaRho = m_originRadius;
115  } else {
116  // originZ = 0.;
117  GlobalPoint myTmp(0., 0., 0.);
118  vertex = myTmp;
119  deltaRho = 1.;
120  deltaZVertex = 15.;
121  }
122 
123  // get the jet direction
125  e.getByToken(token_jet, h_jets);
126 
130  e.getByToken(token_measurementTracker, hmte);
131  measurementTracker = hmte.product();
132  }
133 
134  const auto& bfield = es.getData(token_bfield);
135  const auto& msmaker = es.getData(token_msmaker);
136 
137  for (const reco::Candidate& myJet : *h_jets) {
138  GlobalVector jetVector(myJet.momentum().x(), myJet.momentum().y(), myJet.momentum().z());
139  // GlobalPoint vertex(0, 0, originZ);
140  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(jetVector,
141  vertex,
142  m_ptMin,
143  deltaRho,
144  deltaZVertex,
145  m_deltaEta,
146  m_deltaPhi,
147  bfield,
148  &msmaker,
149  true,
152  m_searchOpt));
153  }
154 
155  return result;
156  }
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
T const * product() const
Definition: Handle.h:70
constexpr bool isUninitialized() const noexcept
Definition: EDGetToken.h:104
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::EDGetTokenT< MeasurementTrackerEvent > token_measurementTracker
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > token_bfield
edm::EDGetTokenT< reco::VertexCollection > token_vertex
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_howToUseMeasurementTracker
edm::ESGetToken< MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord > token_msmaker
edm::EDGetTokenT< reco::CandidateView > token_jet

Member Data Documentation

◆ m_deltaEta

float TauRegionalPixelSeedGenerator::m_deltaEta
private

Definition at line 162 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

◆ m_deltaPhi

float TauRegionalPixelSeedGenerator::m_deltaPhi
private

Definition at line 163 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

◆ m_halfLength

float TauRegionalPixelSeedGenerator::m_halfLength
private

Definition at line 161 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

◆ m_howToUseMeasurementTracker

RectangularEtaPhiTrackingRegion::UseMeasurementTracker TauRegionalPixelSeedGenerator::m_howToUseMeasurementTracker
private

Definition at line 169 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

◆ m_originRadius

float TauRegionalPixelSeedGenerator::m_originRadius
private

Definition at line 160 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

◆ m_ptMin

float TauRegionalPixelSeedGenerator::m_ptMin
private

Definition at line 159 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

◆ m_searchOpt

bool TauRegionalPixelSeedGenerator::m_searchOpt
private

Definition at line 170 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

◆ token_bfield

edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> TauRegionalPixelSeedGenerator::token_bfield
private

Definition at line 167 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions().

◆ token_jet

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

Definition at line 165 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

◆ token_measurementTracker

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

Definition at line 166 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

◆ token_msmaker

edm::ESGetToken<MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord> TauRegionalPixelSeedGenerator::token_msmaker
private

Definition at line 168 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions().

◆ token_vertex

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

Definition at line 164 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().