CMS 3D CMS Logo

List of all members | Public Member Functions | Private Types | Private Attributes | Static Private Attributes
LocalMaximumSeedFinder Class Referencefinal
Inheritance diagram for LocalMaximumSeedFinder:
SeedFinderBase

Public Member Functions

void findSeeds (const edm::Handle< reco::PFRecHitCollection > &input, const std::vector< bool > &mask, std::vector< bool > &seedable) override
 
 LocalMaximumSeedFinder (const edm::ParameterSet &conf)
 
 LocalMaximumSeedFinder (const LocalMaximumSeedFinder &)=delete
 
LocalMaximumSeedFinderoperator= (const LocalMaximumSeedFinder &)=delete
 
- Public Member Functions inherited from SeedFinderBase
const std::string & name () const
 
SeedFinderBaseoperator= (const SeedFinderBase &)=delete
 
 SeedFinderBase (const edm::ParameterSet &conf)
 
 SeedFinderBase (const SeedFinderBase &)=delete
 
virtual ~SeedFinderBase ()=default
 

Private Types

typedef std::tuple< std::vector< int >, std::vector< double >, std::vector< double > > I3tuple
 

Private Attributes

const std::unordered_map< std::string, int > _layerMap
 
const int _nNeighbours
 
std::array< I3tuple, 35 > _thresholds
 

Static Private Attributes

static constexpr double detacut = 0.01
 
static constexpr double dphicut = 0.01
 
static constexpr int layerOffset = 15
 

Detailed Description

Definition at line 12 of file LocalMaximumSeedFinder.cc.

Member Typedef Documentation

◆ I3tuple

typedef std::tuple<std::vector<int>, std::vector<double>, std::vector<double> > LocalMaximumSeedFinder::I3tuple
private

Definition at line 27 of file LocalMaximumSeedFinder.cc.

Constructor & Destructor Documentation

◆ LocalMaximumSeedFinder() [1/2]

LocalMaximumSeedFinder::LocalMaximumSeedFinder ( const edm::ParameterSet conf)

Definition at line 42 of file LocalMaximumSeedFinder.cc.

References PFLayer::ECAL_BARREL, PFLayer::ECAL_ENDCAP, PFLayer::HCAL_BARREL1, PFLayer::HCAL_BARREL2, PFLayer::HCAL_ENDCAP, PFLayer::HF_EM, PFLayer::HF_HAD, createfilelist::int, PFLayer::NONE, PFLayer::PS1, and PFLayer::PS2.

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 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::array< I3tuple, 35 > _thresholds
static constexpr int layerOffset
VParameterSet const & getParameterSetVector(std::string const &name) const
const std::unordered_map< std::string, int > _layerMap
SeedFinderBase(const edm::ParameterSet &conf)

◆ LocalMaximumSeedFinder() [2/2]

LocalMaximumSeedFinder::LocalMaximumSeedFinder ( const LocalMaximumSeedFinder )
delete

Member Function Documentation

◆ findSeeds()

void LocalMaximumSeedFinder::findSeeds ( const edm::Handle< reco::PFRecHitCollection > &  input,
const std::vector< bool > &  mask,
std::vector< bool > &  seedable 
)
overridevirtual

Implements SeedFinderBase.

Definition at line 94 of file LocalMaximumSeedFinder.cc.

References _nNeighbours, _thresholds, funct::abs(), fileCollector::cmp, submitPVResolutionJobs::count, declareDynArray, SiPixelRawToDigiRegional_cfi::deltaPhi, hcalRecHitTable_cff::depth, detacut, dphicut, Exception, PFLayer::HCAL_BARREL1, PFLayer::HCAL_BARREL2, PFLayer::HCAL_ENDCAP, PFLayer::HF_EM, PFLayer::HF_HAD, mps_fire::i, heavyIonCSV_trainingSettings::idx, initDynArray, input, createfilelist::int, dqmiolumiharvest::j, layerOffset, LogDebug, gpuClustering::pixelStatus::mask, eostools::move(), nhits, particleFlowZeroSuppressionECAL_cff::thresholds, and unInitDynArray.

96  {
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 }
#define initDynArray(T, n, x, i)
Definition: DynArray.h:94
constexpr uint32_t mask
Definition: gpuClustering.h:26
static std::string const input
Definition: EdmProvDump.cc:50
std::array< I3tuple, 35 > _thresholds
#define unInitDynArray(T, n, x)
Definition: DynArray.h:88
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static constexpr int layerOffset
static constexpr double detacut
static constexpr double dphicut
#define declareDynArray(T, n, x)
Definition: DynArray.h:91
def move(src, dest)
Definition: eostools.py:511
#define LogDebug(id)

◆ operator=()

LocalMaximumSeedFinder& LocalMaximumSeedFinder::operator= ( const LocalMaximumSeedFinder )
delete

Member Data Documentation

◆ _layerMap

const std::unordered_map<std::string, int> LocalMaximumSeedFinder::_layerMap
private

Definition at line 25 of file LocalMaximumSeedFinder.cc.

◆ _nNeighbours

const int LocalMaximumSeedFinder::_nNeighbours
private

Definition at line 23 of file LocalMaximumSeedFinder.cc.

Referenced by findSeeds().

◆ _thresholds

std::array<I3tuple, 35> LocalMaximumSeedFinder::_thresholds
private

Definition at line 29 of file LocalMaximumSeedFinder.cc.

Referenced by findSeeds().

◆ detacut

constexpr double LocalMaximumSeedFinder::detacut = 0.01
staticprivate

Definition at line 32 of file LocalMaximumSeedFinder.cc.

Referenced by findSeeds().

◆ dphicut

constexpr double LocalMaximumSeedFinder::dphicut = 0.01
staticprivate

Definition at line 33 of file LocalMaximumSeedFinder.cc.

Referenced by findSeeds().

◆ layerOffset

constexpr int LocalMaximumSeedFinder::layerOffset = 15
staticprivate

Definition at line 30 of file LocalMaximumSeedFinder.cc.

Referenced by findSeeds().