CMS 3D CMS Logo

LocalMaximumSeedFinder.cc
Go to the documentation of this file.
5 
6 #include <algorithm>
7 #include <cfloat>
8 #include <tuple>
9 #include <unordered_map>
10 #include <queue>
11 
12 class LocalMaximumSeedFinder final : public SeedFinderBase {
13 public:
17 
19  const std::vector<bool>& mask,
20  std::vector<bool>& seedable) override;
21 
22 private:
23  const int _nNeighbours;
24 
25  const std::unordered_map<std::string, int> _layerMap;
26 
27  typedef std::tuple<std::vector<int>, std::vector<double>, std::vector<double> > I3tuple;
28 
29  std::array<I3tuple, 35> _thresholds;
30  static constexpr int layerOffset = 15;
31 
32  static constexpr double detacut = 0.01;
33  static constexpr double dphicut = 0.01;
34 };
35 
37 
38 namespace {
39  const reco::PFRecHit::Neighbours _noNeighbours(nullptr, 0);
40 }
41 
43  : SeedFinderBase(conf),
44  _nNeighbours(conf.getParameter<int>("nNeighbours")),
45  _layerMap({{"PS2", (int)PFLayer::PS2},
46  {"PS1", (int)PFLayer::PS1},
47  {"ECAL_ENDCAP", (int)PFLayer::ECAL_ENDCAP},
48  {"ECAL_BARREL", (int)PFLayer::ECAL_BARREL},
49  {"NONE", (int)PFLayer::NONE},
50  {"HCAL_BARREL1", (int)PFLayer::HCAL_BARREL1},
51  {"HCAL_BARREL2_RING0", (int)PFLayer::HCAL_BARREL2},
52  // hack to deal with ring1 in HO
53  {"HCAL_BARREL2_RING1", 19},
54  {"HCAL_ENDCAP", (int)PFLayer::HCAL_ENDCAP},
55  {"HF_EM", (int)PFLayer::HF_EM},
56  {"HF_HAD", (int)PFLayer::HF_HAD}}) {
57  const std::vector<edm::ParameterSet>& thresholds = conf.getParameterSetVector("thresholdsByDetector");
58  for (const auto& pset : thresholds) {
59  const std::string& det = pset.getParameter<std::string>("detector");
60 
61  std::vector<int> depths;
62  std::vector<double> thresh_E;
63  std::vector<double> thresh_pT;
64  std::vector<double> thresh_pT2;
65 
66  if (det == std::string("HCAL_BARREL1") || det == std::string("HCAL_ENDCAP")) {
67  depths = pset.getParameter<std::vector<int> >("depths");
68  thresh_E = pset.getParameter<std::vector<double> >("seedingThreshold");
69  thresh_pT = pset.getParameter<std::vector<double> >("seedingThresholdPt");
70  if (thresh_E.size() != depths.size() || thresh_pT.size() != depths.size()) {
71  throw cms::Exception("InvalidGatheringThreshold") << "gatheringThresholds mismatch with the numbers of depths";
72  }
73  } else {
74  depths.push_back(0);
75  thresh_E.push_back(pset.getParameter<double>("seedingThreshold"));
76  thresh_pT.push_back(pset.getParameter<double>("seedingThresholdPt"));
77  }
78 
79  for (unsigned int i = 0; i < thresh_pT.size(); ++i) {
80  thresh_pT2.push_back(thresh_pT[i] * thresh_pT[i]);
81  }
82 
83  auto entry = _layerMap.find(det);
84  if (entry == _layerMap.end()) {
85  throw cms::Exception("InvalidDetectorLayer") << "Detector layer : " << det << " is not in the list of recognized"
86  << " detector layers!";
87  }
88 
89  _thresholds[entry->second + layerOffset] = std::make_tuple(depths, thresh_E, thresh_pT2);
90  }
91 }
92 
93 // the starting state of seedable is all false!
95  const std::vector<bool>& mask,
96  std::vector<bool>& seedable) {
97  auto nhits = input->size();
98  initDynArray(bool, nhits, usable, true);
99  //need to run over energy sorted rechits
100  declareDynArray(float, nhits, energies);
101  unInitDynArray(int, nhits, qst); // queue storage
102  auto cmp = [&](int i, int j) { return energies[i] < energies[j]; };
103  std::priority_queue<int, DynArray<int>, decltype(cmp)> ordered_hits(cmp, std::move(qst));
104 
105  for (unsigned i = 0; i < nhits; ++i) {
106  if (!mask[i])
107  continue; // cannot seed masked objects
108  auto const& maybeseed = (*input)[i];
109  energies[i] = maybeseed.energy();
110  int seedlayer = (int)maybeseed.layer();
111  if (seedlayer == PFLayer::HCAL_BARREL2 && std::abs(maybeseed.positionREP().eta()) > 0.34) {
112  seedlayer = 19;
113  }
114  auto const& thresholds = _thresholds[seedlayer + layerOffset];
115 
116  double thresholdE = 0.;
117  double thresholdPT2 = 0.;
118 
119  for (unsigned int j = 0; j < (std::get<2>(thresholds)).size(); ++j) {
120  int depth = std::get<0>(thresholds)[j];
121  if ((seedlayer == PFLayer::HCAL_BARREL1 && maybeseed.depth() == depth) ||
122  (seedlayer == PFLayer::HCAL_ENDCAP && maybeseed.depth() == depth) ||
123  (seedlayer != PFLayer::HCAL_BARREL1 && seedlayer != PFLayer::HCAL_ENDCAP)) {
124  thresholdE = std::get<1>(thresholds)[j];
125  thresholdPT2 = std::get<2>(thresholds)[j];
126  }
127  }
128 
129  if (maybeseed.energy() < thresholdE || maybeseed.pt2() < thresholdPT2)
130  usable[i] = false;
131  if (!usable[i])
132  continue;
133  ordered_hits.push(i);
134  }
135 
136  while (!ordered_hits.empty()) {
137  auto idx = ordered_hits.top();
138  ordered_hits.pop();
139  if (!usable[idx])
140  continue;
141  //get the neighbours of this seed
142  auto const& maybeseed = (*input)[idx];
143  reco::PFRecHit::Neighbours myNeighbours;
144  switch (_nNeighbours) {
145  case -1:
146  myNeighbours = maybeseed.neighbours();
147  break;
148  case 0: // for HF clustering
149  myNeighbours = _noNeighbours;
150  break;
151  case 4:
152  myNeighbours = maybeseed.neighbours4();
153  break;
154  case 8:
155  myNeighbours = maybeseed.neighbours8();
156  break;
157  default:
158  throw cms::Exception("InvalidConfiguration") << "LocalMaximumSeedFinder only accepts nNeighbors = {-1,0,4,8}";
159  }
160  seedable[idx] = true;
161  for (auto neighbour : myNeighbours) {
162  if (!mask[neighbour])
163  continue;
164  if (energies[neighbour] > energies[idx]) {
165  // std::cout << "how this can be?" << std::endl;
166  seedable[idx] = false;
167  break;
168  }
169  }
170  if (seedable[idx]) {
171  for (auto neighbour : myNeighbours) {
172  //
173  // For HCAL,
174  // even if channel a is a neighbor of channel b, channel b may not be a neighbor of channel a.
175  // So, perform additional checks to ensure making a hit unusable for seeding is safe.
176  int seedlayer = (int)maybeseed.layer();
177  switch (seedlayer) {
180  case PFLayer::HF_EM: // with the current HF setting, we won't see this case
181  // but this can come in if we change _nNeighbours for HF.
182  case PFLayer::HF_HAD: // same as above
183  // HO has only one depth and eta-phi segmentation is regular, so no need to make this check
184  auto const& nei = (*input)[neighbour];
185  if (maybeseed.depth() != nei.depth())
186  continue; // masking is done only if the neighbor is on the same depth layer as the seed
187  if (std::abs(deltaPhi(maybeseed.positionREP().phi(), nei.positionREP().phi())) > dphicut &&
188  std::abs(maybeseed.positionREP().eta() - nei.positionREP().eta()) > detacut)
189  continue; // masking is done only if the neighbor is on the swiss-cross w.r.t. the seed
190  break;
191  }
192 
193  usable[neighbour] = false;
194 
195  } // for-loop
196  }
197  }
198 
199  LogDebug("LocalMaximumSeedFinder") << " found " << std::count(seedable.begin(), seedable.end(), true) << " seeds";
200 }
SeedFinderBase
Definition: SeedFinderBase.h:9
mps_fire.i
i
Definition: mps_fire.py:428
input
static const std::string input
Definition: EdmProvDump.cc:48
MessageLogger.h
LocalMaximumSeedFinder::layerOffset
static constexpr int layerOffset
Definition: LocalMaximumSeedFinder.cc:30
SeedFinderBase.h
particleFlowZeroSuppressionECAL_cff.thresholds
thresholds
Definition: particleFlowZeroSuppressionECAL_cff.py:31
LocalMaximumSeedFinder::detacut
static constexpr double detacut
Definition: LocalMaximumSeedFinder.cc:32
deltaPhi.h
mps_splice.entry
entry
Definition: mps_splice.py:68
PFLayer::HCAL_ENDCAP
Definition: PFLayer.h:37
LocalMaximumSeedFinder::dphicut
static constexpr double dphicut
Definition: LocalMaximumSeedFinder.cc:33
initDynArray
#define initDynArray(T, n, x, i)
Definition: DynArray.h:94
edm::Handle< reco::PFRecHitCollection >
PFLayer::ECAL_BARREL
Definition: PFLayer.h:33
heavyIonCSV_trainingSettings.idx
idx
Definition: heavyIonCSV_trainingSettings.py:5
PFLayer::PS1
Definition: PFLayer.h:31
PFLayer::HCAL_BARREL2
Definition: PFLayer.h:36
PFLayer::HF_EM
Definition: PFLayer.h:38
SiPixelRawToDigiRegional_cfi.deltaPhi
deltaPhi
Definition: SiPixelRawToDigiRegional_cfi.py:9
LocalMaximumSeedFinder::_layerMap
const std::unordered_map< std::string, int > _layerMap
Definition: LocalMaximumSeedFinder.cc:25
PFLayer::HCAL_BARREL1
Definition: PFLayer.h:35
declareDynArray
#define declareDynArray(T, n, x)
Definition: DynArray.h:91
submitPVResolutionJobs.count
count
Definition: submitPVResolutionJobs.py:352
PFLayer::NONE
Definition: PFLayer.h:34
nhits
Definition: HIMultiTrackSelector.h:42
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition: PluginFactory.h:124
LEDCalibrationChannels.depth
depth
Definition: LEDCalibrationChannels.py:65
PFLayer::HF_HAD
Definition: PFLayer.h:39
LocalMaximumSeedFinder
Definition: LocalMaximumSeedFinder.cc:12
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:233
edm::ParameterSet
Definition: ParameterSet.h:47
edmplugin::PluginFactory
Definition: PluginFactory.h:34
LocalMaximumSeedFinder::_nNeighbours
const int _nNeighbours
Definition: LocalMaximumSeedFinder.cc:23
createfilelist.int
int
Definition: createfilelist.py:10
LocalMaximumSeedFinder::operator=
LocalMaximumSeedFinder & operator=(const LocalMaximumSeedFinder &)=delete
LocalMaximumSeedFinder::findSeeds
void findSeeds(const edm::Handle< reco::PFRecHitCollection > &input, const std::vector< bool > &mask, std::vector< bool > &seedable) override
Definition: LocalMaximumSeedFinder.cc:94
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
DynArray.h
reco::PFRecHit::Neighbours
Definition: PFRecHit.h:39
HLT_FULL_cff.depths
depths
Definition: HLT_FULL_cff.py:13426
eostools.move
def move(src, dest)
Definition: eostools.py:511
fileCollector.cmp
cmp
Definition: fileCollector.py:125
LocalMaximumSeedFinder::LocalMaximumSeedFinder
LocalMaximumSeedFinder(const edm::ParameterSet &conf)
Definition: LocalMaximumSeedFinder.cc:42
Exception
Definition: hltDiff.cc:245
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
LocalMaximumSeedFinder::I3tuple
std::tuple< std::vector< int >, std::vector< double >, std::vector< double > > I3tuple
Definition: LocalMaximumSeedFinder.cc:27
PFLayer::ECAL_ENDCAP
Definition: PFLayer.h:32
PFLayer::PS2
Definition: PFLayer.h:30
unInitDynArray
#define unInitDynArray(T, n, x)
Definition: DynArray.h:88
LocalMaximumSeedFinder::_thresholds
std::array< I3tuple, 35 > _thresholds
Definition: LocalMaximumSeedFinder.cc:29
muonDTDigis_cfi.pset
pset
Definition: muonDTDigis_cfi.py:27
findQualityFiles.size
size
Write out results.
Definition: findQualityFiles.py:443