CMS 3D CMS Logo

HGCalClusteringImpl.cc
Go to the documentation of this file.
1 #include <unordered_set>
2 #include <unordered_map>
6 
7 //class constructor
9  : siliconSeedThreshold_(conf.getParameter<double>("seeding_threshold_silicon")),
10  siliconTriggerCellThreshold_(conf.getParameter<double>("clustering_threshold_silicon")),
11  scintillatorSeedThreshold_(conf.getParameter<double>("seeding_threshold_scintillator")),
12  scintillatorTriggerCellThreshold_(conf.getParameter<double>("clustering_threshold_scintillator")),
13  dr_(conf.getParameter<double>("dR_cluster")),
14  clusteringAlgorithmType_(conf.getParameter<string>("clusterType")),
15  calibSF_(conf.getParameter<double>("calibSF_cluster")),
16  layerWeights_(conf.getParameter<std::vector<double>>("layerWeights")),
17  applyLayerWeights_(conf.getParameter<bool>("applyLayerCalibration")) {
18  edm::LogInfo("HGCalClusterParameters") << "C2d Clustering Algorithm selected : " << clusteringAlgorithmType_;
19  edm::LogInfo("HGCalClusterParameters") << "C2d silicon seeding Thr: " << siliconSeedThreshold_;
20  edm::LogInfo("HGCalClusterParameters") << "C2d silicon clustering Thr: " << siliconTriggerCellThreshold_;
21  edm::LogInfo("HGCalClusterParameters") << "C2d scintillator seeding Thr: " << scintillatorSeedThreshold_;
22  edm::LogInfo("HGCalClusterParameters") << "C2d scintillator clustering Thr: " << scintillatorTriggerCellThreshold_;
23  edm::LogInfo("HGCalClusterParameters") << "C2d global calibration factor: " << calibSF_;
24 }
25 
26 /* dR-algorithms */
28  const l1t::HGCalCluster& clu,
29  double distXY) const {
30  DetId tcDetId(tc.detId());
31  DetId cluDetId(clu.detId());
32  if ((triggerTools_.layer(tcDetId) != triggerTools_.layer(cluDetId)) || (tcDetId.subdetId() != cluDetId.subdetId()) ||
33  (triggerTools_.zside(tcDetId) != triggerTools_.zside(cluDetId))) {
34  return false;
35  }
36  if (clu.distance((tc)) < distXY) {
37  return true;
38  }
39  return false;
40 }
41 
44  bool isSeed[triggerCellsPtrs.size()];
45 
46  /* search for cluster seeds */
47  int itc(0);
48  for (std::vector<edm::Ptr<l1t::HGCalTriggerCell>>::const_iterator tc = triggerCellsPtrs.begin();
49  tc != triggerCellsPtrs.end();
50  ++tc, ++itc) {
51  double seedThreshold = ((*tc)->subdetId() == HGCHEB ? scintillatorSeedThreshold_ : siliconSeedThreshold_);
52  isSeed[itc] = ((*tc)->mipPt() > seedThreshold) ? true : false;
53  }
54 
55  /* clustering the TCs */
56  std::vector<l1t::HGCalCluster> clustersTmp;
57 
58  itc = 0;
59  for (std::vector<edm::Ptr<l1t::HGCalTriggerCell>>::const_iterator tc = triggerCellsPtrs.begin();
60  tc != triggerCellsPtrs.end();
61  ++tc, ++itc) {
63  if ((*tc)->mipPt() < threshold) {
64  continue;
65  }
66 
67  /*
68  searching for TC near the center of the cluster
69  ToBeFixed: if a tc is not a seed, but could be potencially be part of a cluster generated by a late seed,
70  the tc will not be clusterized
71  */
72 
73  double minDist = dr_;
74  int targetClu = -1;
75 
76  for (unsigned iclu = 0; iclu < clustersTmp.size(); iclu++) {
77  if (!this->isPertinent(**tc, clustersTmp.at(iclu), dr_))
78  continue;
79 
80  double d = clustersTmp.at(iclu).distance(**tc);
81  if (d < minDist) {
82  minDist = d;
83  targetClu = int(iclu);
84  }
85  }
86 
87  if (targetClu < 0 && isSeed[itc])
88  clustersTmp.emplace_back(*tc);
89  else if (targetClu >= 0)
90  clustersTmp.at(targetClu).addConstituent(*tc);
91  }
92 
93  /* store clusters in the persistent collection */
94  clusters.resize(0, clustersTmp.size());
95  for (unsigned i(0); i < clustersTmp.size(); ++i) {
96  calibratePt(clustersTmp.at(i));
97  clusters.set(0, i, clustersTmp.at(i));
98  }
99 }
100 
101 /* NN-algorithms */
102 
103 /* storing trigger cells into vector per layer and per endcap */
105  const std::vector<edm::Ptr<l1t::HGCalTriggerCell>>& triggerCellsPtrs,
106  std::array<std::vector<std::vector<edm::Ptr<l1t::HGCalTriggerCell>>>, kNSides_>& reshuffledTriggerCells) {
107  for (const auto& tc : triggerCellsPtrs) {
108  DetId tcDetId(tc->detId());
109  int endcap = triggerTools_.zside(tcDetId) == -1 ? 0 : 1;
110  unsigned layer = triggerTools_.layerWithOffset(tc->detId());
111 
112  reshuffledTriggerCells[endcap][layer - 1].emplace_back(tc);
113  }
114 }
115 
116 /* merge clusters that have common neighbors */
118  const l1t::HGCalCluster& secondary_cluster) const {
119  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& pertinentTC = secondary_cluster.constituents();
120 
121  for (const auto& id_tc : pertinentTC) {
122  main_cluster.addConstituent(id_tc.second);
123  }
124 }
125 
128  const HGCalTriggerGeometryBase& triggerGeometry) {
129  /* declaring the clusters vector */
130  std::vector<l1t::HGCalCluster> clustersTmp;
131 
132  // map TC id -> cluster index in clustersTmp
133  std::unordered_map<uint32_t, unsigned> cluNNmap;
134 
135  /* loop over the trigger-cells */
136  for (const auto& tc_ptr : reshuffledTriggerCells) {
137  double threshold =
139  if (tc_ptr->mipPt() < threshold) {
140  continue;
141  }
142 
143  // Check if the neighbors of that TC are already included in a cluster
144  // If this is the case, add the TC to the first (arbitrary) neighbor cluster
145  // Otherwise create a new cluster
146  bool createNewC2d(true);
147  const auto neighbors = triggerGeometry.getNeighborsFromTriggerCell(tc_ptr->detId());
148  for (const auto neighbor : neighbors) {
149  auto tc_cluster_itr = cluNNmap.find(neighbor);
150  if (tc_cluster_itr != cluNNmap.end()) {
151  createNewC2d = false;
152  if (tc_cluster_itr->second < clustersTmp.size()) {
153  clustersTmp.at(tc_cluster_itr->second).addConstituent(tc_ptr);
154  // map TC id to the existing cluster
155  cluNNmap.emplace(tc_ptr->detId(), tc_cluster_itr->second);
156  } else {
157  throw cms::Exception("HGCTriggerUnexpected")
158  << "Trying to access a non-existing cluster. But it should exist...\n";
159  }
160  break;
161  }
162  }
163  if (createNewC2d) {
164  clustersTmp.emplace_back(tc_ptr);
165  clustersTmp.back().setValid(true);
166  // map TC id to the cluster index (size - 1)
167  cluNNmap.emplace(tc_ptr->detId(), clustersTmp.size() - 1);
168  }
169  }
170 
171  /* declaring the vector with possible clusters merged */
172  // Merge neighbor clusters together
173  for (auto& cluster1 : clustersTmp) {
174  // If the cluster has been merged into another one, skip it
175  if (!cluster1.valid())
176  continue;
177  // Fill a set containing all TC included in the clusters
178  // as well as all neighbor TC
179  std::unordered_set<uint32_t> cluTcSet;
180  for (const auto& tc_clu1 : cluster1.constituents()) {
181  cluTcSet.insert(tc_clu1.second->detId());
182  const auto neighbors = triggerGeometry.getNeighborsFromTriggerCell(tc_clu1.second->detId());
183  for (const auto neighbor : neighbors) {
184  cluTcSet.insert(neighbor);
185  }
186  }
187 
188  for (auto& cluster2 : clustersTmp) {
189  // If the cluster has been merged into another one, skip it
190  if (cluster1.detId() == cluster2.detId())
191  continue;
192  if (!cluster2.valid())
193  continue;
194  // Check if the TC in clu2 are in clu1 or its neighbors
195  // If yes, merge the second cluster into the first one
196  for (const auto& tc_clu2 : cluster2.constituents()) {
197  if (cluTcSet.find(tc_clu2.second->detId()) != cluTcSet.end()) {
198  mergeClusters(cluster1, cluster2);
199  cluTcSet.insert(tc_clu2.second->detId());
200  const auto neighbors = triggerGeometry.getNeighborsFromTriggerCell(tc_clu2.second->detId());
201  for (const auto neighbor : neighbors) {
202  cluTcSet.insert(neighbor);
203  }
204  cluster2.setValid(false);
205  break;
206  }
207  }
208  }
209  }
210 
211  /* store clusters in the persistent collection */
212  // only if the cluster contain a TC above the seed threshold
213  for (auto& cluster : clustersTmp) {
214  if (!cluster.valid())
215  continue;
216  bool saveInCollection(false);
217  for (const auto& id_tc : cluster.constituents()) {
218  /* threshold in transverse-mip */
219  double seedThreshold = (id_tc.second->subdetId() == HGCHEB ? scintillatorSeedThreshold_ : siliconSeedThreshold_);
220  if (id_tc.second->mipPt() > seedThreshold) {
221  saveInCollection = true;
222  break;
223  }
224  }
225  if (saveInCollection) {
226  calibratePt(cluster);
227  clusters.push_back(0, cluster);
228  }
229  }
230 }
231 
234  const HGCalTriggerGeometryBase& triggerGeometry) {
235  std::array<std::vector<std::vector<edm::Ptr<l1t::HGCalTriggerCell>>>, kNSides_> reshuffledTriggerCells;
237  for (unsigned side = 0; side < kNSides_; side++) {
238  reshuffledTriggerCells[side].resize(layers);
239  }
240  triggerCellReshuffling(triggerCellsPtrs, reshuffledTriggerCells);
241 
242  for (unsigned iec = 0; iec < kNSides_; ++iec) {
243  for (unsigned il = 0; il < layers; ++il) {
244  NNKernel(reshuffledTriggerCells[iec][il], clusters, triggerGeometry);
245  }
246  }
247 }
248 
249 /*** FW-algorithms ***/
252  const HGCalTriggerGeometryBase& triggerGeometry) {
253  bool isSeed[triggerCellsPtrs.size()];
254  std::vector<unsigned> seedPositions;
255  seedPositions.reserve(triggerCellsPtrs.size());
256 
257  /* search for cluster seeds */
258  int itc(0);
259  for (std::vector<edm::Ptr<l1t::HGCalTriggerCell>>::const_iterator tc = triggerCellsPtrs.begin();
260  tc != triggerCellsPtrs.end();
261  ++tc, ++itc) {
262  double seedThreshold = ((*tc)->subdetId() == HGCHEB ? scintillatorSeedThreshold_ : siliconSeedThreshold_);
263 
264  /* decide if is a seed, if yes store the position into of triggerCellsPtrs */
265  isSeed[itc] = ((*tc)->mipPt() > seedThreshold) ? true : false;
266  if (isSeed[itc]) {
267  seedPositions.push_back(itc);
268 
269  /* remove tc from the seed vector if is a NN of an other seed*/
270  for (auto pos : seedPositions) {
271  if (((*tc)->position() - triggerCellsPtrs[pos]->position()).mag() < dr_) {
272  if (this->areTCneighbour((*tc)->detId(), triggerCellsPtrs[pos]->detId(), triggerGeometry)) {
273  isSeed[itc] = false;
274  seedPositions.pop_back();
275  }
276  }
277  }
278  }
279  }
280 
281  /* clustering the TCs */
282  std::vector<l1t::HGCalCluster> clustersTmp;
283 
284  // every seed generates a cluster
285  clustersTmp.reserve(seedPositions.size());
286  for (auto pos : seedPositions) {
287  clustersTmp.emplace_back(triggerCellsPtrs[pos]);
288  }
289 
290  /* add the tc to the clusters */
291  itc = 0;
292  for (std::vector<edm::Ptr<l1t::HGCalTriggerCell>>::const_iterator tc = triggerCellsPtrs.begin();
293  tc != triggerCellsPtrs.end();
294  ++tc, ++itc) {
295  /* get the correct threshold for the different part of the detector */
297 
298  /* continue if not passing the threshold */
299  if ((*tc)->mipPt() < threshold)
300  continue;
301  if (isSeed[itc])
302  continue; //No sharing of seeds between clusters (TBC)
303 
304  /* searching for TC near the centre of the cluster */
305  std::vector<unsigned> tcPertinentClusters;
306  unsigned iclu(0);
307 
308  for (const auto& clu : clustersTmp) {
309  if (this->isPertinent(**tc, clu, dr_)) {
310  tcPertinentClusters.push_back(iclu);
311  }
312  ++iclu;
313  }
314 
315  if (tcPertinentClusters.empty()) {
316  continue;
317  } else if (tcPertinentClusters.size() == 1) {
318  clustersTmp.at(tcPertinentClusters.at(0)).addConstituent(*tc);
319  } else {
320  /* calculate the fractions */
321  double totMipt = 0;
322  for (auto clu : tcPertinentClusters) {
323  totMipt += clustersTmp.at(clu).seedMipPt();
324  }
325 
326  for (auto clu : tcPertinentClusters) {
327  double seedMipt = clustersTmp.at(clu).seedMipPt();
328  clustersTmp.at(clu).addConstituent(*tc, true, seedMipt / totMipt);
329  }
330  }
331  }
332 
333  /* store clusters in the persistent collection */
334  clusters.resize(0, clustersTmp.size());
335  for (unsigned i(0); i < clustersTmp.size(); ++i) {
336  this->removeUnconnectedTCinCluster(clustersTmp.at(i), triggerGeometry);
337  calibratePt(clustersTmp.at(i));
338  clusters.set(0, i, clustersTmp.at(i));
339  }
340 }
341 
343  uint32_t detIDb,
344  const HGCalTriggerGeometryBase& triggerGeometry) {
345  const auto neighbors = triggerGeometry.getNeighborsFromTriggerCell(detIDa);
346 
347  if (neighbors.find(detIDb) != neighbors.end())
348  return true;
349 
350  return false;
351 }
352 
354  const HGCalTriggerGeometryBase& triggerGeometry) {
355  /* get the constituents and the centre of the seed tc (considered as the first of the constituents) */
356  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& constituents = cluster.constituents();
357  Basic3DVector<float> seedCentre(constituents.at(cluster.detId())->position());
358 
359  /* distances from the seed */
360  vector<pair<edm::Ptr<l1t::HGCalTriggerCell>, float>> distances;
361  for (const auto& id_tc : constituents) {
362  Basic3DVector<float> tcCentre(id_tc.second->position());
363  float distance = (seedCentre - tcCentre).mag();
364  distances.push_back(pair<edm::Ptr<l1t::HGCalTriggerCell>, float>(id_tc.second, distance));
365  }
366 
367  /* sorting (needed in order to be sure that we are skipping any tc) */
368  /* FIXME: better sorting needed!!! */
369  std::sort(distances.begin(), distances.end(), distanceSorter);
370 
371  /* checking if the tc is connected to the seed */
372  bool toRemove[constituents.size()];
373  toRemove[0] = false; // this is the seed
374  for (unsigned itc = 1; itc < distances.size(); itc++) {
375  /* get the tc under study */
376  toRemove[itc] = true;
377  const edm::Ptr<l1t::HGCalTriggerCell>& tcToStudy = distances[itc].first;
378 
379  /* compare with the tc in the cluster */
380  for (unsigned itc_ref = 0; itc_ref < itc; itc_ref++) {
381  if (!toRemove[itc_ref]) {
382  if (areTCneighbour(tcToStudy->detId(), distances.at(itc_ref).first->detId(), triggerGeometry)) {
383  toRemove[itc] = false;
384  break;
385  }
386  }
387  }
388  }
389 
390  /* remove the unconnected TCs */
391  for (unsigned i = 0; i < distances.size(); i++) {
392  if (toRemove[i])
393  cluster.removeConstituent(distances.at(i).first);
394  }
395 }
396 
398  double calibPt = 0.;
399 
400  if (applyLayerWeights_) {
401  unsigned layerN = triggerTools_.layerWithOffset(cluster.detId());
402 
403  if (layerWeights_.at(layerN) == 0.) {
404  throw cms::Exception("BadConfiguration")
405  << "2D cluster energy forced to 0 by calibration coefficients.\n"
406  << "The configuration should be changed. "
407  << "Discarded layers should be defined in hgcalTriggerGeometryESProducer.TriggerGeometry.DisconnectedLayers "
408  "and not with calibration coefficients = 0\n";
409  }
410 
411  calibPt = layerWeights_.at(layerN) * cluster.mipPt();
412 
413  }
414 
415  else {
416  calibPt = cluster.pt() * calibSF_;
417  }
418 
419  math::PtEtaPhiMLorentzVector calibP4(calibPt, cluster.eta(), cluster.phi(), 0.);
420 
421  cluster.setP4(calibP4);
422 }
void addConstituent(const edm::Ptr< C > &c, bool updateCentre=true, float fraction=1.)
Definition: HGCalClusterT.h:54
double pt() const final
transverse momentum
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:49
bool distanceSorter(pair< edm::Ptr< l1t::HGCalTriggerCell >, float > i, pair< edm::Ptr< l1t::HGCalTriggerCell >, float > j)
void calibratePt(l1t::HGCalCluster &cluster)
uint32_t detId() const
Definition: HGCalClusterT.h:92
virtual geom_set getNeighborsFromTriggerCell(const unsigned trigger_cell_det_id) const =0
std::vector< double > layerWeights_
bool isPertinent(const l1t::HGCalTriggerCell &tc, const l1t::HGCalCluster &clu, double distXY) const
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
int zside(const DetId &) const
void clusterizeDRNN(const std::vector< edm::Ptr< l1t::HGCalTriggerCell >> &triggerCellsPtrs, l1t::HGCalClusterBxCollection &clusters, const HGCalTriggerGeometryBase &triggerGeometry)
void removeConstituent(const edm::Ptr< C > &c, bool updateCentre=true)
Definition: HGCalClusterT.h:72
static constexpr unsigned kNSides_
uint32_t detId() const
void clusterizeNN(const std::vector< edm::Ptr< l1t::HGCalTriggerCell >> &triggerCellsPtrs, l1t::HGCalClusterBxCollection &clusters, const HGCalTriggerGeometryBase &triggerGeometry)
std::string clusteringAlgorithmType_
void removeUnconnectedTCinCluster(l1t::HGCalCluster &cluster, const HGCalTriggerGeometryBase &triggerGeometry)
void triggerCellReshuffling(const std::vector< edm::Ptr< l1t::HGCalTriggerCell >> &triggerCellsPtrs, std::array< std::vector< std::vector< edm::Ptr< l1t::HGCalTriggerCell >>>, kNSides_ > &reshuffledTriggerCells)
d
Definition: ztail.py:151
unsigned layerWithOffset(const DetId &) const
double distance(const l1t::HGCalTriggerCell &tc) const
Definition: HGCalClusterT.h:97
Log< level::Info, false > LogInfo
Definition: DetId.h:17
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
HGCalTriggerTools triggerTools_
unsigned layers(ForwardSubdetector type) const
HGCalClusteringImpl(const edm::ParameterSet &conf)
void NNKernel(const std::vector< edm::Ptr< l1t::HGCalTriggerCell >> &reshuffledTriggerCells, l1t::HGCalClusterBxCollection &clusters, const HGCalTriggerGeometryBase &triggerGeometry)
void clusterizeDR(const std::vector< edm::Ptr< l1t::HGCalTriggerCell >> &triggerCellsPtrs, l1t::HGCalClusterBxCollection &clusters)
unsigned layer(const DetId &) const
bool areTCneighbour(uint32_t detIDa, uint32_t detIDb, const HGCalTriggerGeometryBase &triggerGeometry)
static int position[264][3]
Definition: ReadPGInfo.cc:289
double scintillatorTriggerCellThreshold_
double phi() const final
momentum azimuthal angle
double mipPt() const
Definition: HGCalClusterT.h:90
void setP4(const LorentzVector &p4) final
set 4-momentum
void mergeClusters(l1t::HGCalCluster &main_cluster, const l1t::HGCalCluster &secondary_cluster) const
double eta() const final
momentum pseudorapidity