|
|
Go to the documentation of this file. 1 #ifndef PointSeededTrackingRegionsProducer_h
2 #define PointSeededTrackingRegionsProducer_h
22 #include <TLorentzVector.h>
57 if (modeString ==
"BeamSpotFixed")
59 else if (modeString ==
"BeamSpotSigma")
61 else if (modeString ==
"VerticesFixed")
63 else if (modeString ==
"VerticesSigma")
66 edm::LogError(
"PointSeededTrackingRegionsProducer") <<
"Unknown mode string: " << modeString;
97 if (regPSet.
exists(
"searchOpt"))
110 <<
"nSigmaZBeamSpot must be positive for BeamSpotSigma mode!";
120 descPoints.
add<std::vector<double>>(
"eta", {0.});
121 descPoints.
add<std::vector<double>>(
"phi", {0.});
125 desc.add<
int>(
"maxNRegions", 10);
128 desc.add<
int>(
"maxNVertices", 1);
130 desc.add<
double>(
"ptMin", 0.9);
131 desc.add<
double>(
"originRadius", 0.2);
132 desc.add<
double>(
"zErrorBeamSpot", 24.2);
133 desc.add<
double>(
"deltaEta", 0.5);
134 desc.add<
double>(
"deltaPhi", 0.5);
135 desc.add<
bool>(
"precise",
true);
137 desc.add<
double>(
"nSigmaZVertex", 3.);
138 desc.add<
double>(
"zErrorVetex", 0.2);
139 desc.add<
double>(
"nSigmaZBeamSpot", 4.);
144 desc.add<
bool>(
"searchOpt",
false);
152 descriptions.
add(
"pointSeededTrackingRegion", descRegion);
156 std::vector<std::unique_ptr<TrackingRegion>>
result;
168 std::vector<std::pair<GlobalPoint, float>> origins;
172 origins.push_back(std::make_pair(
178 for (reco::VertexCollection::const_iterator iv =
vertices->begin(),
ev =
vertices->end();
181 if (iv->isFake() || !iv->isValid())
184 origins.push_back(std::make_pair(
GlobalPoint(iv->x(), iv->y(), iv->z()),
189 if (origins.empty()) {
190 origins.push_back(std::make_pair(
215 result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(direction,
229 edm::LogInfo(
"PointSeededTrackingRegionsProducer") <<
"produced " << n_regions <<
" regions";
ParameterDescriptionBase * add(U const &iLabel, T const &value)
T const * product() const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::vector< std::unique_ptr< TrackingRegion > > regions(const edm::Event &e, const edm::EventSetup &es) const override
std::vector< Vertex > VertexCollection
collection of Vertex objects
constexpr bool isUninitialized() const noexcept
U second(std::pair< T, U > const &p)
Log< level::Info, false > LogInfo
Sin< T >::type sin(const T &t)
std::vector< double > etaPoints
Cos< T >::type cos(const T &t)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Geom::Theta< T > theta() const
Global3DPoint GlobalPoint
~PointSeededTrackingRegionsProducer() override
bool exists(std::string const ¶meterName) const
checks if a parameter exists
std::vector< double > phiPoints
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot
Tan< T >::type tan(const T &t)
Log< level::Error, false > LogError
PointSeededTrackingRegionsProducer(const edm::ParameterSet &conf, edm::ConsumesCollector &&iC)
edm::EDGetTokenT< MeasurementTrackerEvent > token_measurementTracker
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_whereToUseMeasurementTracker
T getParameter(std::string const &) const
static UseMeasurementTracker stringToUseMeasurementTracker(const std::string &name)
edm::EDGetTokenT< reco::VertexCollection > token_vertex