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 {
16  edm::LogInfo("HGCalClusterParameters") << "C2d Clustering Algorithm selected : " << clusteringAlgorithmType_ ;
17  edm::LogInfo("HGCalClusterParameters") << "C2d silicon seeding Thr: " << siliconSeedThreshold_ ;
18  edm::LogInfo("HGCalClusterParameters") << "C2d silicon clustering Thr: " << siliconTriggerCellThreshold_ ;
19  edm::LogInfo("HGCalClusterParameters") << "C2d scintillator seeding Thr: " << scintillatorSeedThreshold_ ;
20  edm::LogInfo("HGCalClusterParameters") << "C2d scintillator clustering Thr: " << scintillatorTriggerCellThreshold_ ;
21 
22 }
23 
24 
25 /* dR-algorithms */
27  const l1t::HGCalCluster & clu,
28  double distXY ) const
29 {
30 
31  HGCalDetId tcDetId( tc.detId() );
32  HGCalDetId cluDetId( clu.detId() );
33  if( (tcDetId.layer() != cluDetId.layer()) ||
34  (tcDetId.subdetId() != cluDetId.subdetId()) ||
35  (tcDetId.zside() != cluDetId.zside()) ){
36  return false;
37  }
38  if ( clu.distance((tc)) < distXY ){
39  return true;
40  }
41  return false;
42 
43 }
44 
45 
48  ){
49 
50  bool isSeed[triggerCellsPtrs.size()];
51 
52  /* search for cluster seeds */
53  int itc(0);
54  for( edm::PtrVector<l1t::HGCalTriggerCell>::const_iterator tc = triggerCellsPtrs.begin(); tc != triggerCellsPtrs.end(); ++tc,++itc ){
55  double seedThreshold = ((*tc)->subdetId()==HGCHEB ? scintillatorSeedThreshold_ : siliconSeedThreshold_);
56  isSeed[itc] = ( (*tc)->mipPt() > seedThreshold) ? true : false;
57  }
58 
59  /* clustering the TCs */
60  std::vector<l1t::HGCalCluster> clustersTmp;
61 
62  itc=0;
63  for( edm::PtrVector<l1t::HGCalTriggerCell>::const_iterator tc = triggerCellsPtrs.begin(); tc != triggerCellsPtrs.end(); ++tc,++itc ){
65  if( (*tc)->mipPt() < threshold ){
66  continue;
67  }
68 
69  /* searching for TC near the center of the cluster */
70  int iclu=0;
71  vector<int> tcPertinentClusters;
72  for( const auto& clu : clustersTmp){
73  if( this->isPertinent(**tc, clu, dr_) ){
74  tcPertinentClusters.push_back(iclu);
75  }
76  ++iclu;
77  }
78  if( tcPertinentClusters.empty() && isSeed[itc] ){
79  clustersTmp.emplace_back( *tc );
80  }
81  else if ( !tcPertinentClusters.empty() ){
82 
83  unsigned minDist(300);
84  unsigned targetClu(0);
85 
86  for( int iclu : tcPertinentClusters){
87  double d = clustersTmp.at(iclu).distance(**tc);
88  if( d < minDist ){
89  minDist = d;
90  targetClu = iclu;
91  }
92  }
93 
94  clustersTmp.at(targetClu).addConstituent( *tc );
95 
96  }
97  }
98 
99  /* store clusters in the persistent collection */
100  clusters.resize(0, clustersTmp.size());
101  for( unsigned i(0); i<clustersTmp.size(); ++i ){
102  clusters.set( 0, i, clustersTmp.at(i) );
103  }
104 
105 }
106 
107 
108 
109 /* NN-algorithms */
110 
111 /* storing trigger cells into vector per layer and per endcap */
113  std::array< std::array<std::vector<edm::Ptr<l1t::HGCalTriggerCell>>, kLayers_>,kNSides_> & reshuffledTriggerCells
114  ){
115 
116  for( const auto& tc : triggerCellsPtrs ){
117  int endcap = tc->zside() == -1 ? 0 : 1 ;
118  HGCalDetId tcDetId( tc->detId() );
119  int subdet = tcDetId.subdetId();
120  int layer = -1;
121 
122  if( subdet == HGCEE ){
123  layer = tc->layer();
124  }
125  else if( subdet == HGCHEF ){
126  layer = tc->layer() + kLayersEE_;
127  }
128  else if( subdet == HGCHEB ){
129  layer = tc->layer() + kLayersEE_ + kLayersFH_;
130  }
131 
132  reshuffledTriggerCells[endcap][layer-1].emplace_back(tc);
133 
134  }
135 
136 }
137 
138 
139 /* merge clusters that have common neighbors */
141  const l1t::HGCalCluster & secondary_cluster ) const
142 {
143 
144  const edm::PtrVector<l1t::HGCalTriggerCell>& pertinentTC = secondary_cluster.constituents();
145 
146  for( edm::PtrVector<l1t::HGCalTriggerCell>::iterator tc = pertinentTC.begin(); tc != pertinentTC.end(); ++tc ){
147  main_cluster.addConstituent(*tc);
148  }
149 
150 }
151 
152 
153 void HGCalClusteringImpl::NNKernel( const std::vector<edm::Ptr<l1t::HGCalTriggerCell>> & reshuffledTriggerCells,
155  const HGCalTriggerGeometryBase & triggerGeometry
156  ){
157 
158  /* declaring the clusters vector */
159  std::vector<l1t::HGCalCluster> clustersTmp;
160 
161  // map TC id -> cluster index in clustersTmp
162  std::unordered_map<uint32_t, unsigned> cluNNmap;
163 
164  /* loop over the trigger-cells */
165  for( const auto& tc_ptr : reshuffledTriggerCells ){
167  if( tc_ptr->mipPt() < threshold ){
168  continue;
169  }
170 
171  // Check if the neighbors of that TC are already included in a cluster
172  // If this is the case, add the TC to the first (arbitrary) neighbor cluster
173  // Otherwise create a new cluster
174  bool createNewC2d(true);
175  const auto neighbors = triggerGeometry.getNeighborsFromTriggerCell(tc_ptr->detId());
176  for( const auto neighbor : neighbors ){
177  auto tc_cluster_itr = cluNNmap.find(neighbor);
178  if(tc_cluster_itr!=cluNNmap.end()){
179  createNewC2d = false;
180  if( tc_cluster_itr->second < clustersTmp.size()){
181  clustersTmp.at(tc_cluster_itr->second).addConstituent(tc_ptr);
182  // map TC id to the existing cluster
183  cluNNmap.emplace(tc_ptr->detId(), tc_cluster_itr->second);
184  }
185  else{
186  throw cms::Exception("HGCTriggerUnexpected")
187  << "Trying to access a non-existing cluster. But it should exist...\n";
188  }
189  break;
190  }
191  }
192  if(createNewC2d){
193  clustersTmp.emplace_back(tc_ptr);
194  clustersTmp.back().setValid(true);
195  // map TC id to the cluster index (size - 1)
196  cluNNmap.emplace(tc_ptr->detId(), clustersTmp.size()-1);
197  }
198  }
199 
200  /* declaring the vector with possible clusters merged */
201  // Merge neighbor clusters together
202  for( auto& cluster1 : clustersTmp){
203  // If the cluster has been merged into another one, skip it
204  if( !cluster1.valid() ) continue;
205  // Fill a set containing all TC included in the clusters
206  // as well as all neighbor TC
207  std::unordered_set<uint32_t> cluTcSet;
208  for(const auto& tc_clu1 : cluster1.constituents()){
209  cluTcSet.insert( tc_clu1->detId() );
210  const auto neighbors = triggerGeometry.getNeighborsFromTriggerCell( tc_clu1->detId() );
211  for(const auto neighbor : neighbors){
212  cluTcSet.insert( neighbor );
213  }
214  }
215 
216  for( auto& cluster2 : clustersTmp ){
217  // If the cluster has been merged into another one, skip it
218  if( cluster1.detId()==cluster2.detId() ) continue;
219  if( !cluster2.valid() ) continue;
220  // Check if the TC in clu2 are in clu1 or its neighbors
221  // If yes, merge the second cluster into the first one
222  for(const auto& tc_clu2 : cluster2.constituents()){
223  if( cluTcSet.find(tc_clu2->detId())!=cluTcSet.end() ){
224  mergeClusters( cluster1, cluster2 );
225  cluTcSet.insert( tc_clu2->detId() );
226  const auto neighbors = triggerGeometry.getNeighborsFromTriggerCell( tc_clu2->detId() );
227  for(const auto neighbor : neighbors){
228  cluTcSet.insert( neighbor );
229  }
230  cluster2.setValid(false);
231  break;
232  }
233  }
234  }
235  }
236 
237  /* store clusters in the persistent collection */
238  // only if the cluster contain a TC above the seed threshold
239  for( auto& cluster : clustersTmp ){
240  if( !cluster.valid() ) continue;
241  bool saveInCollection(false);
242  for( const auto& tc_ptr : cluster.constituents() ){
243  /* threshold in transverse-mip */
244  double seedThreshold = (tc_ptr->subdetId()==HGCHEB ? scintillatorSeedThreshold_ : siliconSeedThreshold_);
245  if( tc_ptr->mipPt() > seedThreshold ){
246  saveInCollection = true;
247  break;
248  }
249  }
250  if(saveInCollection){
251  clusters.push_back( 0, cluster );
252  }
253  }
254 }
255 
256 
259  const HGCalTriggerGeometryBase & triggerGeometry
260  ){
261 
262  std::array< std::array< std::vector<edm::Ptr<l1t::HGCalTriggerCell> >,kLayers_>,kNSides_> reshuffledTriggerCells;
263  triggerCellReshuffling( triggerCellsPtrs, reshuffledTriggerCells );
264 
265  for(unsigned iec=0; iec<kNSides_; ++iec){
266  for(unsigned il=0; il<kLayers_; ++il){
267  NNKernel( reshuffledTriggerCells[iec][il], clusters, triggerGeometry );
268  }
269  }
270 
271 }
272 
double distance(const l1t::HGCalTriggerCell &tc) const
Definition: HGCalClusterT.h:99
size_type size() const
Size of the RefVector.
Definition: PtrVectorBase.h:74
void clusterizeDR(const edm::PtrVector< l1t::HGCalTriggerCell > &triggerCellsPtrs, l1t::HGCalClusterBxCollection &clusters)
virtual geom_set getNeighborsFromTriggerCell(const unsigned trigger_cell_det_id) const =0
void mergeClusters(l1t::HGCalCluster &main_cluster, const l1t::HGCalCluster &secondary_cluster) const
static const unsigned kLayersFH_
void clusterizeNN(const edm::PtrVector< l1t::HGCalTriggerCell > &triggerCellsPtrs, l1t::HGCalClusterBxCollection &clusters, const HGCalTriggerGeometryBase &triggerGeometry)
const_iterator begin() const
Definition: PtrVector.h:129
uint32_t detId() const
bool isPertinent(const l1t::HGCalTriggerCell &tc, const l1t::HGCalCluster &clu, double distXY) const
void triggerCellReshuffling(const edm::PtrVector< l1t::HGCalTriggerCell > &triggerCellsPtrs, std::array< std::array< std::vector< edm::Ptr< l1t::HGCalTriggerCell >>, kLayers_ >, kNSides_ > &reshuffledTriggerCells)
static const unsigned kLayersEE_
std::string clusteringAlgorithmType_
const_iterator end() const
Definition: PtrVector.h:134
uint32_t detId() const
Definition: HGCalClusterT.h:96
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
void set(int bx, unsigned i, const T &object)
static const unsigned kNSides_
void addConstituent(const edm::Ptr< C > &c)
Definition: HGCalClusterT.h:54
HGCalClusteringImpl(const edm::ParameterSet &conf)
void resize(int bx, unsigned size)
const edm::PtrVector< C > & constituents() const
Definition: HGCalClusterT.h:49
void NNKernel(const std::vector< edm::Ptr< l1t::HGCalTriggerCell >> &reshuffledTriggerCells, l1t::HGCalClusterBxCollection &clusters, const HGCalTriggerGeometryBase &triggerGeometry)
static const unsigned kLayers_
double scintillatorTriggerCellThreshold_
void push_back(int bx, T object)