1 #ifndef __InitialClusteringStepBase_H__ 2 #define __InitialClusteringStepBase_H__ 14 #include <unordered_map> 43 _algoName(conf.getParameter<
std::string>(
"algoName")) {
44 const std::vector<edm::ParameterSet>&
thresholds = conf.getParameterSetVector(
"thresholdsByDetector");
45 for (
const auto&
pset : thresholds) {
49 std::vector<double> thresh_E;
50 std::vector<double> thresh_pT;
51 std::vector<double> thresh_pT2;
54 depths =
pset.getParameter<std::vector<int> >(
"depths");
55 thresh_E =
pset.getParameter<std::vector<double> >(
"gatheringThreshold");
56 thresh_pT =
pset.getParameter<std::vector<double> >(
"gatheringThresholdPt");
57 if (thresh_E.size() != depths.size() || thresh_pT.size() != depths.size()) {
59 <<
"gatheringThresholds mismatch with the numbers of depths";
63 thresh_E.push_back(
pset.getParameter<
double>(
"gatheringThreshold"));
64 thresh_pT.push_back(
pset.getParameter<
double>(
"gatheringThresholdPt"));
67 for (
unsigned int i = 0;
i < thresh_pT.size(); ++
i) {
68 thresh_pT2.push_back(thresh_pT[
i] * thresh_pT[
i]);
71 auto entry = _layerMap.find(det);
72 if (
entry == _layerMap.end()) {
74 <<
"Detector layer : " << det <<
" is not in the list of recognized" 75 <<
" detector layers!";
77 _thresholds.emplace(_layerMap.find(det)->second, std::make_tuple(depths, thresh_E, thresh_pT2));
83 ICSB& operator=(
const ICSB&) =
delete;
90 const std::vector<bool>& mask,
91 const std::vector<bool>& seeds,
95 o <<
"InitialClusteringStep with algo \"" << _algoName <<
"\" located " << _nSeeds <<
" seeds and built " 96 << _nClustersFound <<
" clusters from those seeds. ";
100 void reset() { _nSeeds = _nClustersFound = 0; }
109 typedef std::tuple<std::vector<int>, std::vector<double>, std::vector<double> >
I3tuple;
const std::unordered_map< std::string, int > _layerMap
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
reco::PFRecHitRef makeRefhit(const edm::Handle< reco::PFRecHitCollection > &h, const unsigned i) const
InitialClusteringStepBase ICSB
virtual void updateEvent(const edm::Event &)
const std::string _algoName
edm::Ref< PFRecHitCollection > PFRecHitRef
persistent reference to PFRecHit objects
edmplugin::PluginFactory< InitialClusteringStepBase *(const edm::ParameterSet &, edm::ConsumesCollector &)> InitialClusteringStepFactory
virtual void update(const edm::EventSetup &)
std::tuple< std::vector< int >, std::vector< double >, std::vector< double > > I3tuple
InitialClusteringStepBase(const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
std::ostream & operator<<(std::ostream &ost, const HLTGlobalStatus &hlt)
Formatted printout of trigger tbale.
std::unordered_map< int, I3tuple > _thresholds
std::ostream & operator<<(std::ostream &o) const