CMS 3D CMS Logo

Basic2DGenericTopoClusterizer.cc
Go to the documentation of this file.
4 
7 
8 public:
10  : InitialClusteringStepBase(conf, cc), _useCornerCells(conf.getParameter<bool>("useCornerCells")) {}
11  ~Basic2DGenericTopoClusterizer() override = default;
12  Basic2DGenericTopoClusterizer(const B2DGT&) = delete;
13  B2DGT& operator=(const B2DGT&) = delete;
14 
16  const std::vector<bool>&,
17  const std::vector<bool>&,
18  reco::PFClusterCollection&) override;
19 
20 private:
21  const bool _useCornerCells;
23  const std::vector<bool>&, // masked rechits
24  unsigned int, //present rechit
25  std::vector<bool>&, // hit usage state
26  reco::PFCluster&); // the topocluster
27 };
28 
30 
31 #ifdef PFLOW_DEBUG
32 #define LOGVERB(x) edm::LogVerbatim(x)
33 #define LOGWARN(x) edm::LogWarning(x)
34 #define LOGERR(x) edm::LogError(x)
35 #define LOGDRESSED(x) edm::LogInfo(x)
36 #else
37 #define LOGVERB(x) LogTrace(x)
38 #define LOGWARN(x) edm::LogWarning(x)
39 #define LOGERR(x) edm::LogError(x)
40 #define LOGDRESSED(x) LogDebug(x)
41 #endif
42 
44  const std::vector<bool>& rechitMask,
45  const std::vector<bool>& seedable,
47  auto const& hits = *input;
48  std::vector<bool> used(hits.size(), false);
49  std::vector<unsigned int> seeds;
50 
51  // get the seeds and sort them descending in energy
52  seeds.reserve(hits.size());
53  for (unsigned int i = 0; i < hits.size(); ++i) {
54  if (!rechitMask[i] || !seedable[i] || used[i])
55  continue;
56  seeds.emplace_back(i);
57  }
58  // maxHeap would be better
59  std::sort(
60  seeds.begin(), seeds.end(), [&](unsigned int i, unsigned int j) { return hits[i].energy() > hits[j].energy(); });
61 
63  for (auto seed : seeds) {
64  if (!rechitMask[seed] || !seedable[seed] || used[seed])
65  continue;
66  temp.reset();
67  buildTopoCluster(input, rechitMask, seed, used, temp);
68  if (!temp.recHitFractions().empty())
69  output.push_back(temp);
70  }
71 }
72 
74  const std::vector<bool>& rechitMask,
75  unsigned int kcell,
76  std::vector<bool>& used,
77  reco::PFCluster& topocluster) {
78  auto const& cell = (*input)[kcell];
79  int cell_layer = (int)cell.layer();
80  if (cell_layer == PFLayer::HCAL_BARREL2 && std::abs(cell.positionREP().eta()) > 0.34) {
81  cell_layer *= 100;
82  }
83 
84  auto const& thresholds = _thresholds.find(cell_layer)->second;
85  double thresholdE = 0.;
86  double thresholdPT2 = 0.;
87 
88  for (unsigned int j = 0; j < (std::get<1>(thresholds)).size(); ++j) {
89  int depth = std::get<0>(thresholds)[j];
90 
91  if ((cell_layer == PFLayer::HCAL_BARREL1 && cell.depth() == depth) ||
92  (cell_layer == PFLayer::HCAL_ENDCAP && cell.depth() == depth) ||
93  (cell_layer != PFLayer::HCAL_BARREL1 && cell_layer != PFLayer::HCAL_ENDCAP)) {
94  thresholdE = std::get<1>(thresholds)[j];
95  thresholdPT2 = std::get<2>(thresholds)[j];
96  }
97  }
98 
99  if (cell.energy() < thresholdE || cell.pt2() < thresholdPT2) {
100  LOGDRESSED("GenericTopoCluster::buildTopoCluster()")
101  << "RecHit " << cell.detId() << " with enegy " << cell.energy() << " GeV was rejected!." << std::endl;
102  return;
103  }
104 
105  auto k = kcell;
106  used[k] = true;
107  auto ref = makeRefhit(input, k);
108  topocluster.addRecHitFraction(reco::PFRecHitFraction(ref, 1.0));
109 
110  auto const& neighbours = (_useCornerCells ? cell.neighbours8() : cell.neighbours4());
111 
112  for (auto nb : neighbours) {
113  if (used[nb] || !rechitMask[nb]) {
114  LOGDRESSED("GenericTopoCluster::buildTopoCluster()")
115  << " RecHit " << cell.detId() << "\'s"
116  << " neighbor RecHit " << input->at(nb).detId() << " with enegy " << input->at(nb).energy()
117  << " GeV was rejected!"
118  << " Reasons : " << used[nb] << " (used) " << !rechitMask[nb] << " (masked)." << std::endl;
119  continue;
120  }
121  buildTopoCluster(input, rechitMask, nb, used, topocluster);
122  }
123 }
reco::PFRecHitRef makeRefhit(const edm::Handle< reco::PFRecHitCollection > &h, const unsigned i) const
~Basic2DGenericTopoClusterizer() override=default
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:42
void buildClusters(const edm::Handle< reco::PFRecHitCollection > &, const std::vector< bool > &, const std::vector< bool > &, reco::PFClusterCollection &) override
void buildTopoCluster(const edm::Handle< reco::PFRecHitCollection > &, const std::vector< bool > &, unsigned int, std::vector< bool > &, reco::PFCluster &)
Fraction of a PFRecHit (rechits can be shared between several PFCluster&#39;s)
Basic2DGenericTopoClusterizer B2DGT
static std::string const input
Definition: EdmProvDump.cc:47
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Basic2DGenericTopoClusterizer(const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
void addRecHitFraction(const reco::PFRecHitFraction &frac)
add a given fraction of the rechit
Definition: PFCluster.cc:33
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
#define DEFINE_EDM_PLUGIN(factory, type, name)
B2DGT & operator=(const B2DGT &)=delete
#define LOGDRESSED(x)
std::unordered_map< int, I3tuple > _thresholds