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 37 of file TauRegionalPixelSeedGenerator.h.

Constructor & Destructor Documentation

◆ TauRegionalPixelSeedGenerator()

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

Definition at line 39 of file TauRegionalPixelSeedGenerator.h.

References edm::ParameterSet::exists(), edm::ParameterSet::existsAs(), edm::ParameterSet::getParameter(), HLT_2022v12_cff::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.

41  edm::LogInfo("TauRegionalPixelSeedGenerator") << "Enter the TauRegionalPixelSeedGenerator";
42 
43  edm::ParameterSet regionPSet = conf.getParameter<edm::ParameterSet>("RegionPSet");
44 
45  m_ptMin = regionPSet.getParameter<double>("ptMin");
46  m_originRadius = regionPSet.getParameter<double>("originRadius");
47  m_halfLength = regionPSet.getParameter<double>("originHalfLength");
48  m_deltaEta = regionPSet.getParameter<double>("deltaEtaRegion");
49  m_deltaPhi = regionPSet.getParameter<double>("deltaPhiRegion");
51  token_vertex = iC.consumes<reco::VertexCollection>(regionPSet.getParameter<edm::InputTag>("vertexSrc"));
52  if (regionPSet.exists("searchOpt")) {
53  m_searchOpt = regionPSet.getParameter<bool>("searchOpt");
54  } else {
55  m_searchOpt = false;
56  }
58  regionPSet.getParameter<std::string>("howToUseMeasurementTracker"));
60  // temporary until everything migrated to InputTag
62  if (regionPSet.existsAs<edm::InputTag>("measurementTrackerName")) {
63  tag = regionPSet.getParameter<edm::InputTag>("measurementTrackerName");
64  } else {
65  tag = edm::InputTag(regionPSet.getParameter<std::string>("measurementTrackerName"));
66  }
67 
69  }
70  }
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 72 of file TauRegionalPixelSeedGenerator.h.

72 {}

Member Function Documentation

◆ fillDescriptions()

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

Definition at line 74 of file TauRegionalPixelSeedGenerator.h.

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

74  {
76 
77  desc.add<double>("ptMin", 5.0);
78  desc.add<double>("originHalfLength", 0.2);
79  desc.add<double>("originRadius", 0.2);
80  desc.add<double>("deltaEtaRegion", 0.1);
81  desc.add<double>("deltaPhiRegion", 0.1);
82  desc.add<edm::InputTag>("JetSrc", edm::InputTag("icone5Tau1"));
83  desc.add<edm::InputTag>("vertexSrc", edm::InputTag("pixelVertices"));
84  desc.add<bool>("searchOpt", false);
85 
86  desc.add<std::string>("howToUseMeasurementTracker", "ForSiStrips");
87 
88  // allow both InputTag and string for the moment, use InputTag as the default
90  "measurementTrackerName", edm::InputTag("MeasurementTrackerEvent"), true) xor
91  edm::ParameterDescription<std::string>("measurementTrackerName", "MeasurementTrackerEvent", true));
92 
93  // Only for backwards-compatibility
95  descRegion.add<edm::ParameterSetDescription>("RegionPSet", desc);
96 
97  descriptions.add("tauRegionalPixelSeedTrackingRegions", descRegion);
98  }
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 100 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.

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

Referenced by regions(), and TauRegionalPixelSeedGenerator().

◆ m_deltaPhi

float TauRegionalPixelSeedGenerator::m_deltaPhi
private

Definition at line 164 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

◆ m_halfLength

float TauRegionalPixelSeedGenerator::m_halfLength
private

Definition at line 162 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

◆ m_howToUseMeasurementTracker

RectangularEtaPhiTrackingRegion::UseMeasurementTracker TauRegionalPixelSeedGenerator::m_howToUseMeasurementTracker
private

Definition at line 170 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

◆ m_originRadius

float TauRegionalPixelSeedGenerator::m_originRadius
private

Definition at line 161 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

◆ m_ptMin

float TauRegionalPixelSeedGenerator::m_ptMin
private

Definition at line 160 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

◆ m_searchOpt

bool TauRegionalPixelSeedGenerator::m_searchOpt
private

Definition at line 171 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

◆ token_bfield

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

Definition at line 168 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions().

◆ token_jet

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

Definition at line 166 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

◆ token_measurementTracker

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

Definition at line 167 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().

◆ token_msmaker

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

Definition at line 169 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions().

◆ token_vertex

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

Definition at line 165 of file TauRegionalPixelSeedGenerator.h.

Referenced by regions(), and TauRegionalPixelSeedGenerator().