CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
HGCalClusteringImpl Class Reference

#include <HGCalClusteringImpl.h>

Public Member Functions

void clusterizeDR (const std::vector< edm::Ptr< l1t::HGCalTriggerCell >> &triggerCellsPtrs, l1t::HGCalClusterBxCollection &clusters)
 
void clusterizeDRNN (const std::vector< edm::Ptr< l1t::HGCalTriggerCell >> &triggerCellsPtrs, l1t::HGCalClusterBxCollection &clusters, const HGCalTriggerGeometryBase &triggerGeometry)
 
void clusterizeNN (const std::vector< edm::Ptr< l1t::HGCalTriggerCell >> &triggerCellsPtrs, l1t::HGCalClusterBxCollection &clusters, const HGCalTriggerGeometryBase &triggerGeometry)
 
void eventSetup (const edm::EventSetup &es)
 
 HGCalClusteringImpl (const edm::ParameterSet &conf)
 
bool isPertinent (const l1t::HGCalTriggerCell &tc, const l1t::HGCalCluster &clu, double distXY) const
 
void mergeClusters (l1t::HGCalCluster &main_cluster, const l1t::HGCalCluster &secondary_cluster) const
 
void NNKernel (const std::vector< edm::Ptr< l1t::HGCalTriggerCell >> &reshuffledTriggerCells, l1t::HGCalClusterBxCollection &clusters, const HGCalTriggerGeometryBase &triggerGeometry)
 

Private Member Functions

bool areTCneighbour (uint32_t detIDa, uint32_t detIDb, const HGCalTriggerGeometryBase &triggerGeometry)
 
void calibratePt (l1t::HGCalCluster &cluster)
 
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)
 

Private Attributes

bool applyLayerWeights_
 
double calibSF_
 
std::string clusteringAlgorithmType_
 
double dr_
 
std::vector< double > layerWeights_
 
double scintillatorSeedThreshold_
 
double scintillatorTriggerCellThreshold_
 
double siliconSeedThreshold_
 
double siliconTriggerCellThreshold_
 
HGCalTriggerTools triggerTools_
 

Static Private Attributes

static const unsigned kNSides_ = 2
 

Detailed Description

Definition at line 18 of file HGCalClusteringImpl.h.

Constructor & Destructor Documentation

HGCalClusteringImpl::HGCalClusteringImpl ( const edm::ParameterSet conf)

Definition at line 8 of file HGCalClusteringImpl.cc.

References calibSF_, clusteringAlgorithmType_, scintillatorSeedThreshold_, scintillatorTriggerCellThreshold_, siliconSeedThreshold_, and siliconTriggerCellThreshold_.

8  :
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 {
19  edm::LogInfo("HGCalClusterParameters") << "C2d Clustering Algorithm selected : " << clusteringAlgorithmType_ ;
20  edm::LogInfo("HGCalClusterParameters") << "C2d silicon seeding Thr: " << siliconSeedThreshold_ ;
21  edm::LogInfo("HGCalClusterParameters") << "C2d silicon clustering Thr: " << siliconTriggerCellThreshold_ ;
22  edm::LogInfo("HGCalClusterParameters") << "C2d scintillator seeding Thr: " << scintillatorSeedThreshold_ ;
23  edm::LogInfo("HGCalClusterParameters") << "C2d scintillator clustering Thr: " << scintillatorTriggerCellThreshold_ ;
24  edm::LogInfo("HGCalClusterParameters") << "C2d global calibration factor: " << calibSF_;
25 }
T getParameter(std::string const &) const
std::vector< double > layerWeights_
std::string clusteringAlgorithmType_
double scintillatorTriggerCellThreshold_

Member Function Documentation

bool HGCalClusteringImpl::areTCneighbour ( uint32_t  detIDa,
uint32_t  detIDb,
const HGCalTriggerGeometryBase triggerGeometry 
)
private

Definition at line 378 of file HGCalClusteringImpl.cc.

References HGCalTriggerGeometryBase::getNeighborsFromTriggerCell().

Referenced by clusterizeDRNN(), and removeUnconnectedTCinCluster().

379  {
380 
381  const auto neighbors = triggerGeometry.getNeighborsFromTriggerCell( detIDa );
382 
383  if( neighbors.find( detIDb ) != neighbors.end() ) return true;
384 
385  return false;
386 
387 }
virtual geom_set getNeighborsFromTriggerCell(const unsigned trigger_cell_det_id) const =0
void HGCalClusteringImpl::calibratePt ( l1t::HGCalCluster cluster)
private

Definition at line 441 of file HGCalClusteringImpl.cc.

References applyLayerWeights_, calibSF_, l1t::HGCalClusterT< C >::detId(), reco::LeafCandidate::eta(), Exception, layerWeights_, HGCalTriggerTools::layerWithOffset(), l1t::HGCalClusterT< C >::mipPt(), reco::LeafCandidate::phi(), reco::LeafCandidate::pt(), reco::LeafCandidate::setP4(), and triggerTools_.

Referenced by clusterizeDR(), clusterizeDRNN(), and NNKernel().

441  {
442 
443  double calibPt=0.;
444 
445  if(applyLayerWeights_){
446 
447  unsigned layerN = triggerTools_.layerWithOffset(cluster.detId());
448 
449  if(layerWeights_.at(layerN)==0.){
450  throw cms::Exception("BadConfiguration")
451  <<"2D cluster energy forced to 0 by calibration coefficients.\n"
452  <<"The configuration should be changed. "
453  <<"Discarded layers should be defined in hgcalTriggerGeometryESProducer.TriggerGeometry.DisconnectedLayers and not with calibration coefficients = 0\n";
454  }
455 
456  calibPt = layerWeights_.at(layerN) * cluster.mipPt();
457 
458  }
459 
460  else{
461 
462  calibPt = cluster.pt() * calibSF_;
463 
464  }
465 
466  math::PtEtaPhiMLorentzVector calibP4( calibPt,
467  cluster.eta(),
468  cluster.phi(),
469  0. );
470 
471  cluster.setP4( calibP4 );
472 
473 }
double eta() const final
momentum pseudorapidity
double pt() const final
transverse momentum
unsigned layerWithOffset(const DetId &) const
std::vector< double > layerWeights_
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
double mipPt() const
uint32_t detId() const
HGCalTriggerTools triggerTools_
double phi() const final
momentum azimuthal angle
void setP4(const LorentzVector &p4) final
set 4-momentum
void HGCalClusteringImpl::clusterizeDR ( const std::vector< edm::Ptr< l1t::HGCalTriggerCell >> &  triggerCellsPtrs,
l1t::HGCalClusterBxCollection clusters 
)

Definition at line 49 of file HGCalClusteringImpl.cc.

References calibratePt(), edmIntegrityCheck::d, dr_, HGCHEB, mps_fire::i, isPertinent(), BXVector< T >::resize(), scintillatorSeedThreshold_, scintillatorTriggerCellThreshold_, UEAnalysisJets_cfi::seedThreshold, BXVector< T >::set(), siliconSeedThreshold_, siliconTriggerCellThreshold_, and electronIdCutBased_cfi::threshold.

Referenced by eventSetup().

51  {
52 
53  bool isSeed[triggerCellsPtrs.size()];
54 
55  /* search for cluster seeds */
56  int itc(0);
57  for( std::vector<edm::Ptr<l1t::HGCalTriggerCell>>::const_iterator tc = triggerCellsPtrs.begin(); tc != triggerCellsPtrs.end(); ++tc,++itc ){
58  double seedThreshold = ((*tc)->subdetId()==HGCHEB ? scintillatorSeedThreshold_ : siliconSeedThreshold_);
59  isSeed[itc] = ( (*tc)->mipPt() > seedThreshold) ? true : false;
60  }
61 
62  /* clustering the TCs */
63  std::vector<l1t::HGCalCluster> clustersTmp;
64 
65  itc=0;
66  for( std::vector<edm::Ptr<l1t::HGCalTriggerCell>>::const_iterator tc = triggerCellsPtrs.begin(); tc != triggerCellsPtrs.end(); ++tc,++itc ){
68  if( (*tc)->mipPt() < threshold ){
69  continue;
70  }
71 
72  /*
73  searching for TC near the center of the cluster
74  ToBeFixed: if a tc is not a seed, but could be potencially be part of a cluster generated by a late seed,
75  the tc will not be clusterized
76  */
77  int iclu=0;
78  vector<int> tcPertinentClusters;
79  for( const auto& clu : clustersTmp){
80  if( this->isPertinent(**tc, clu, dr_) ){
81  tcPertinentClusters.push_back(iclu);
82  }
83  ++iclu;
84  }
85  if( tcPertinentClusters.empty() && isSeed[itc] ){
86  clustersTmp.emplace_back( *tc );
87  }
88  else if ( !tcPertinentClusters.empty() ){
89 
90  unsigned minDist(300);
91  unsigned targetClu(0);
92 
93  for( int iclu : tcPertinentClusters){
94  double d = clustersTmp.at(iclu).distance(**tc);
95  if( d < minDist ){
96  minDist = d;
97  targetClu = iclu;
98  }
99  }
100 
101  clustersTmp.at(targetClu).addConstituent( *tc );
102 
103  }
104  }
105 
106  /* store clusters in the persistent collection */
107  clusters.resize(0, clustersTmp.size());
108  for( unsigned i(0); i<clustersTmp.size(); ++i ){
109  calibratePt(clustersTmp.at(i));
110  clusters.set( 0, i, clustersTmp.at(i) );
111  }
112 
113 
114 
115 }
void calibratePt(l1t::HGCalCluster &cluster)
bool isPertinent(const l1t::HGCalTriggerCell &tc, const l1t::HGCalCluster &clu, double distXY) const
void set(int bx, unsigned i, const T &object)
void resize(int bx, unsigned size)
double scintillatorTriggerCellThreshold_
void HGCalClusteringImpl::clusterizeDRNN ( const std::vector< edm::Ptr< l1t::HGCalTriggerCell >> &  triggerCellsPtrs,
l1t::HGCalClusterBxCollection clusters,
const HGCalTriggerGeometryBase triggerGeometry 
)

Definition at line 281 of file HGCalClusteringImpl.cc.

References areTCneighbour(), calibratePt(), dr_, HGCHEB, mps_fire::i, isPertinent(), mag(), removeUnconnectedTCinCluster(), BXVector< T >::resize(), scintillatorSeedThreshold_, scintillatorTriggerCellThreshold_, UEAnalysisJets_cfi::seedThreshold, BXVector< T >::set(), siliconSeedThreshold_, siliconTriggerCellThreshold_, and electronIdCutBased_cfi::threshold.

Referenced by eventSetup().

284  {
285 
286  bool isSeed[triggerCellsPtrs.size()];
287  std::vector<unsigned> seedPositions;
288  seedPositions.reserve( triggerCellsPtrs.size() );
289 
290  /* search for cluster seeds */
291  int itc(0);
292  for( std::vector<edm::Ptr<l1t::HGCalTriggerCell>>::const_iterator tc = triggerCellsPtrs.begin(); tc != triggerCellsPtrs.end(); ++tc,++itc ){
293 
294  double seedThreshold = ((*tc)->subdetId()==HGCHEB ? scintillatorSeedThreshold_ : siliconSeedThreshold_);
295 
296  /* decide if is a seed, if yes store the position into of triggerCellsPtrs */
297  isSeed[itc] = ( (*tc)->mipPt() > seedThreshold) ? true : false;
298  if( isSeed[itc] ) {
299 
300  seedPositions.push_back( itc );
301 
302  /* remove tc from the seed vector if is a NN of an other seed*/
303  for( auto pos : seedPositions ){
304  if( ( (*tc)->position() - triggerCellsPtrs[pos]->position() ).mag() < dr_ ){
305  if( this->areTCneighbour( (*tc)->detId(), triggerCellsPtrs[pos]->detId(), triggerGeometry ) )
306  {
307  isSeed[itc] = false;
308  seedPositions.pop_back();
309  }
310  }
311  }
312  }
313 
314  }
315 
316  /* clustering the TCs */
317  std::vector<l1t::HGCalCluster> clustersTmp;
318 
319  // every seed generates a cluster
320  for( auto pos : seedPositions ) {
321  clustersTmp.emplace_back( triggerCellsPtrs[pos] );
322  }
323 
324  /* add the tc to the clusters */
325  itc=0;
326  for( std::vector<edm::Ptr<l1t::HGCalTriggerCell>>::const_iterator tc = triggerCellsPtrs.begin(); tc != triggerCellsPtrs.end(); ++tc,++itc ){
327 
328  /* get the correct threshold for the different part of the detector */
330 
331  /* continue if not passing the threshold */
332  if( (*tc)->mipPt() < threshold ) continue;
333  if( isSeed[itc] ) continue; //No sharing of seeds between clusters (TBC)
334 
335  /* searching for TC near the centre of the cluster */
336  std::vector<unsigned> tcPertinentClusters;
337  unsigned iclu(0);
338 
339  for ( const auto& clu : clustersTmp ) {
340  if( this->isPertinent(**tc, clu, dr_) ){
341  tcPertinentClusters.push_back( iclu );
342  }
343  ++iclu;
344  }
345 
346  if ( tcPertinentClusters.empty() ) {
347  continue;
348  }
349  else if( tcPertinentClusters.size() == 1 ) {
350  clustersTmp.at( tcPertinentClusters.at(0) ).addConstituent( *tc );
351  }
352  else {
353 
354  /* calculate the fractions */
355  double totMipt = 0;
356  for( auto clu : tcPertinentClusters ){
357  totMipt += clustersTmp.at( clu ).seedMipPt();
358  }
359 
360  for( auto clu : tcPertinentClusters ){
361  double seedMipt = clustersTmp.at( clu ).seedMipPt();
362  clustersTmp.at( clu ).addConstituent( *tc, true, seedMipt/totMipt );
363  }
364  }
365  }
366 
367  /* store clusters in the persistent collection */
368  clusters.resize(0, clustersTmp.size());
369  for( unsigned i(0); i<clustersTmp.size(); ++i ){
370  this->removeUnconnectedTCinCluster( clustersTmp.at(i), triggerGeometry );
371  calibratePt( clustersTmp.at(i) );
372  clusters.set( 0, i, clustersTmp.at(i) );
373  }
374 
375 }
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
void calibratePt(l1t::HGCalCluster &cluster)
bool isPertinent(const l1t::HGCalTriggerCell &tc, const l1t::HGCalCluster &clu, double distXY) const
void removeUnconnectedTCinCluster(l1t::HGCalCluster &cluster, const HGCalTriggerGeometryBase &triggerGeometry)
void set(int bx, unsigned i, const T &object)
void resize(int bx, unsigned size)
bool areTCneighbour(uint32_t detIDa, uint32_t detIDb, const HGCalTriggerGeometryBase &triggerGeometry)
double scintillatorTriggerCellThreshold_
void HGCalClusteringImpl::clusterizeNN ( const std::vector< edm::Ptr< l1t::HGCalTriggerCell >> &  triggerCellsPtrs,
l1t::HGCalClusterBxCollection clusters,
const HGCalTriggerGeometryBase triggerGeometry 
)

Definition at line 257 of file HGCalClusteringImpl.cc.

References ForwardEmpty, kNSides_, LayerTriplets::layers(), HGCalTriggerTools::layers(), NNKernel(), triggerCellReshuffling(), and triggerTools_.

Referenced by eventSetup().

260  {
261 
262  std::array< std::vector< std::vector<edm::Ptr<l1t::HGCalTriggerCell>>>,kNSides_> reshuffledTriggerCells;
264  for(unsigned side=0; side<kNSides_; side++)
265  {
266  reshuffledTriggerCells[side].resize(layers);
267  }
268  triggerCellReshuffling( triggerCellsPtrs, reshuffledTriggerCells );
269 
270  for(unsigned iec=0; iec<kNSides_; ++iec){
271  for(unsigned il=0; il<layers; ++il){
272  NNKernel( reshuffledTriggerCells[iec][il], clusters, triggerGeometry );
273  }
274  }
275 
276 }
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
void triggerCellReshuffling(const std::vector< edm::Ptr< l1t::HGCalTriggerCell >> &triggerCellsPtrs, std::array< std::vector< std::vector< edm::Ptr< l1t::HGCalTriggerCell >>>, kNSides_ > &reshuffledTriggerCells)
static const unsigned kNSides_
unsigned layers(ForwardSubdetector type) const
HGCalTriggerTools triggerTools_
void NNKernel(const std::vector< edm::Ptr< l1t::HGCalTriggerCell >> &reshuffledTriggerCells, l1t::HGCalClusterBxCollection &clusters, const HGCalTriggerGeometryBase &triggerGeometry)
void HGCalClusteringImpl::eventSetup ( const edm::EventSetup es)
inline
bool HGCalClusteringImpl::isPertinent ( const l1t::HGCalTriggerCell tc,
const l1t::HGCalCluster clu,
double  distXY 
) const

Definition at line 29 of file HGCalClusteringImpl.cc.

References l1t::HGCalTriggerCell::detId(), l1t::HGCalClusterT< C >::detId(), and l1t::HGCalClusterT< C >::distance().

Referenced by clusterizeDR(), clusterizeDRNN(), and eventSetup().

32 {
33 
34  HGCalDetId tcDetId( tc.detId() );
35  HGCalDetId cluDetId( clu.detId() );
36  if( (tcDetId.layer() != cluDetId.layer()) ||
37  (tcDetId.subdetId() != cluDetId.subdetId()) ||
38  (tcDetId.zside() != cluDetId.zside()) ){
39  return false;
40  }
41  if ( clu.distance((tc)) < distXY ){
42  return true;
43  }
44  return false;
45 
46 }
double distance(const l1t::HGCalTriggerCell &tc) const
uint32_t detId() const
uint32_t detId() const
void HGCalClusteringImpl::mergeClusters ( l1t::HGCalCluster main_cluster,
const l1t::HGCalCluster secondary_cluster 
) const

Definition at line 139 of file HGCalClusteringImpl.cc.

References l1t::HGCalClusterT< C >::addConstituent(), and l1t::HGCalClusterT< C >::constituents().

Referenced by eventSetup(), and NNKernel().

141 {
142 
143  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& pertinentTC = secondary_cluster.constituents();
144 
145  for(const auto& id_tc : pertinentTC){
146  main_cluster.addConstituent(id_tc.second);
147  }
148 
149 }
void addConstituent(const edm::Ptr< C > &c, bool updateCentre=true, float fraction=1.)
Definition: HGCalClusterT.h:58
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:52
void HGCalClusteringImpl::NNKernel ( const std::vector< edm::Ptr< l1t::HGCalTriggerCell >> &  reshuffledTriggerCells,
l1t::HGCalClusterBxCollection clusters,
const HGCalTriggerGeometryBase triggerGeometry 
)

Definition at line 152 of file HGCalClusteringImpl.cc.

References calibratePt(), Exception, HGCalTriggerGeometryBase::getNeighborsFromTriggerCell(), HGCHEB, mergeClusters(), BXVector< T >::push_back(), scintillatorSeedThreshold_, scintillatorTriggerCellThreshold_, UEAnalysisJets_cfi::seedThreshold, siliconSeedThreshold_, siliconTriggerCellThreshold_, and electronIdCutBased_cfi::threshold.

Referenced by clusterizeNN(), and eventSetup().

155  {
156 
157  /* declaring the clusters vector */
158  std::vector<l1t::HGCalCluster> clustersTmp;
159 
160  // map TC id -> cluster index in clustersTmp
161  std::unordered_map<uint32_t, unsigned> cluNNmap;
162 
163  /* loop over the trigger-cells */
164  for( const auto& tc_ptr : reshuffledTriggerCells ){
166  if( tc_ptr->mipPt() < threshold ){
167  continue;
168  }
169 
170  // Check if the neighbors of that TC are already included in a cluster
171  // If this is the case, add the TC to the first (arbitrary) neighbor cluster
172  // Otherwise create a new cluster
173  bool createNewC2d(true);
174  const auto neighbors = triggerGeometry.getNeighborsFromTriggerCell(tc_ptr->detId());
175  for( const auto neighbor : neighbors ){
176  auto tc_cluster_itr = cluNNmap.find(neighbor);
177  if(tc_cluster_itr!=cluNNmap.end()){
178  createNewC2d = false;
179  if( tc_cluster_itr->second < clustersTmp.size()){
180  clustersTmp.at(tc_cluster_itr->second).addConstituent(tc_ptr);
181  // map TC id to the existing cluster
182  cluNNmap.emplace(tc_ptr->detId(), tc_cluster_itr->second);
183  }
184  else{
185  throw cms::Exception("HGCTriggerUnexpected")
186  << "Trying to access a non-existing cluster. But it should exist...\n";
187  }
188  break;
189  }
190  }
191  if(createNewC2d){
192  clustersTmp.emplace_back(tc_ptr);
193  clustersTmp.back().setValid(true);
194  // map TC id to the cluster index (size - 1)
195  cluNNmap.emplace(tc_ptr->detId(), clustersTmp.size()-1);
196  }
197  }
198 
199  /* declaring the vector with possible clusters merged */
200  // Merge neighbor clusters together
201  for( auto& cluster1 : clustersTmp){
202  // If the cluster has been merged into another one, skip it
203  if( !cluster1.valid() ) continue;
204  // Fill a set containing all TC included in the clusters
205  // as well as all neighbor TC
206  std::unordered_set<uint32_t> cluTcSet;
207  for(const auto& tc_clu1 : cluster1.constituents()){
208  cluTcSet.insert( tc_clu1.second->detId() );
209  const auto neighbors = triggerGeometry.getNeighborsFromTriggerCell( tc_clu1.second->detId() );
210  for(const auto neighbor : neighbors){
211  cluTcSet.insert( neighbor );
212  }
213  }
214 
215  for( auto& cluster2 : clustersTmp ){
216  // If the cluster has been merged into another one, skip it
217  if( cluster1.detId()==cluster2.detId() ) continue;
218  if( !cluster2.valid() ) continue;
219  // Check if the TC in clu2 are in clu1 or its neighbors
220  // If yes, merge the second cluster into the first one
221  for(const auto& tc_clu2 : cluster2.constituents()){
222  if( cluTcSet.find(tc_clu2.second->detId())!=cluTcSet.end() ){
223  mergeClusters( cluster1, cluster2 );
224  cluTcSet.insert( tc_clu2.second->detId() );
225  const auto neighbors = triggerGeometry.getNeighborsFromTriggerCell( tc_clu2.second->detId() );
226  for(const auto neighbor : neighbors){
227  cluTcSet.insert( neighbor );
228  }
229  cluster2.setValid(false);
230  break;
231  }
232  }
233  }
234  }
235 
236  /* store clusters in the persistent collection */
237  // only if the cluster contain a TC above the seed threshold
238  for( auto& cluster : clustersTmp ){
239  if( !cluster.valid() ) continue;
240  bool saveInCollection(false);
241  for( const auto& id_tc : cluster.constituents() ){
242  /* threshold in transverse-mip */
243  double seedThreshold = (id_tc.second->subdetId()==HGCHEB ? scintillatorSeedThreshold_ : siliconSeedThreshold_);
244  if( id_tc.second->mipPt() > seedThreshold ){
245  saveInCollection = true;
246  break;
247  }
248  }
249  if(saveInCollection){
250  calibratePt(cluster);
251  clusters.push_back( 0, cluster );
252  }
253  }
254 }
virtual geom_set getNeighborsFromTriggerCell(const unsigned trigger_cell_det_id) const =0
void calibratePt(l1t::HGCalCluster &cluster)
void mergeClusters(l1t::HGCalCluster &main_cluster, const l1t::HGCalCluster &secondary_cluster) const
double scintillatorTriggerCellThreshold_
void push_back(int bx, T object)
void HGCalClusteringImpl::removeUnconnectedTCinCluster ( l1t::HGCalCluster cluster,
const HGCalTriggerGeometryBase triggerGeometry 
)
private

Definition at line 390 of file HGCalClusteringImpl.cc.

References areTCneighbour(), l1t::HGCalClusterT< C >::constituents(), l1t::HGCalTriggerCell::detId(), l1t::HGCalClusterT< C >::detId(), SoftLeptonByDistance_cfi::distance, distanceSorter(), mps_fire::i, mag(), position, and l1t::HGCalClusterT< C >::removeConstituent().

Referenced by clusterizeDRNN().

390  {
391 
392  /* get the constituents and the centre of the seed tc (considered as the first of the constituents) */
393  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& constituents = cluster.constituents();
394  Basic3DVector<float> seedCentre( constituents.at(cluster.detId())->position() );
395 
396  /* distances from the seed */
397  vector<pair<edm::Ptr<l1t::HGCalTriggerCell>,float>> distances;
398  for( const auto & id_tc : constituents )
399  {
400  Basic3DVector<float> tcCentre( id_tc.second->position() );
401  float distance = ( seedCentre - tcCentre ).mag();
402  distances.push_back( pair<edm::Ptr<l1t::HGCalTriggerCell>,float>( id_tc.second, distance ) );
403  }
404 
405  /* sorting (needed in order to be sure that we are skipping any tc) */
406  /* FIXME: better sorting needed!!! */
407  std::sort( distances.begin(), distances.end(), distanceSorter );
408 
409 
410  /* checking if the tc is connected to the seed */
411  bool toRemove[constituents.size()];
412  toRemove[0] = false; // this is the seed
413  for( unsigned itc=1; itc<distances.size(); itc++ ){
414 
415  /* get the tc under study */
416  toRemove[itc] = true;
417  const edm::Ptr<l1t::HGCalTriggerCell>& tcToStudy = distances[itc].first;
418 
419  /* compare with the tc in the cluster */
420  for( unsigned itc_ref=0; itc_ref<itc; itc_ref++ ){
421  if( !toRemove[itc_ref] ) {
422  if( areTCneighbour( tcToStudy->detId(), distances.at( itc_ref ).first->detId(), triggerGeometry ) ) {
423  toRemove[itc] = false;
424  break;
425  }
426  }
427  }
428 
429  }
430 
431 
432  /* remove the unconnected TCs */
433  for( unsigned i=0; i<distances.size(); i++){
434  if( toRemove[i] ) cluster.removeConstituent( distances.at( i ).first );
435  }
436 
437 }
bool distanceSorter(pair< edm::Ptr< l1t::HGCalTriggerCell >, float > i, pair< edm::Ptr< l1t::HGCalTriggerCell >, float > j)
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
uint32_t detId() const
void removeConstituent(const edm::Ptr< C > &c, bool updateCentre=true)
Definition: HGCalClusterT.h:82
uint32_t detId() const
bool areTCneighbour(uint32_t detIDa, uint32_t detIDb, const HGCalTriggerGeometryBase &triggerGeometry)
static int position[264][3]
Definition: ReadPGInfo.cc:509
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
Definition: HGCalClusterT.h:52
void HGCalClusteringImpl::triggerCellReshuffling ( const std::vector< edm::Ptr< l1t::HGCalTriggerCell >> &  triggerCellsPtrs,
std::array< std::vector< std::vector< edm::Ptr< l1t::HGCalTriggerCell >>>, kNSides_ > &  reshuffledTriggerCells 
)
private

Definition at line 122 of file HGCalClusteringImpl.cc.

References makeMuonMisalignmentScenario::endcap, HGCalTriggerTools::layerWithOffset(), and triggerTools_.

Referenced by clusterizeNN().

124  {
125 
126  for( const auto& tc : triggerCellsPtrs ){
127  int endcap = tc->zside() == -1 ? 0 : 1 ;
128  HGCalDetId tcDetId( tc->detId() );
129  unsigned layer = triggerTools_.layerWithOffset(tc->detId());
130 
131  reshuffledTriggerCells[endcap][layer-1].emplace_back(tc);
132 
133  }
134 
135 }
unsigned layerWithOffset(const DetId &) const
HGCalTriggerTools triggerTools_

Member Data Documentation

bool HGCalClusteringImpl::applyLayerWeights_
private

Definition at line 70 of file HGCalClusteringImpl.h.

Referenced by calibratePt().

double HGCalClusteringImpl::calibSF_
private

Definition at line 68 of file HGCalClusteringImpl.h.

Referenced by calibratePt(), and HGCalClusteringImpl().

std::string HGCalClusteringImpl::clusteringAlgorithmType_
private

Definition at line 67 of file HGCalClusteringImpl.h.

Referenced by HGCalClusteringImpl().

double HGCalClusteringImpl::dr_
private

Definition at line 66 of file HGCalClusteringImpl.h.

Referenced by clusterizeDR(), and clusterizeDRNN().

const unsigned HGCalClusteringImpl::kNSides_ = 2
staticprivate

Definition at line 21 of file HGCalClusteringImpl.h.

Referenced by clusterizeNN().

std::vector<double> HGCalClusteringImpl::layerWeights_
private

Definition at line 69 of file HGCalClusteringImpl.h.

Referenced by calibratePt().

double HGCalClusteringImpl::scintillatorSeedThreshold_
private

Definition at line 64 of file HGCalClusteringImpl.h.

Referenced by clusterizeDR(), clusterizeDRNN(), HGCalClusteringImpl(), and NNKernel().

double HGCalClusteringImpl::scintillatorTriggerCellThreshold_
private

Definition at line 65 of file HGCalClusteringImpl.h.

Referenced by clusterizeDR(), clusterizeDRNN(), HGCalClusteringImpl(), and NNKernel().

double HGCalClusteringImpl::siliconSeedThreshold_
private

Definition at line 62 of file HGCalClusteringImpl.h.

Referenced by clusterizeDR(), clusterizeDRNN(), HGCalClusteringImpl(), and NNKernel().

double HGCalClusteringImpl::siliconTriggerCellThreshold_
private

Definition at line 63 of file HGCalClusteringImpl.h.

Referenced by clusterizeDR(), clusterizeDRNN(), HGCalClusteringImpl(), and NNKernel().

HGCalTriggerTools HGCalClusteringImpl::triggerTools_
private