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 369 of file HGCalClusteringImpl.cc.

References HGCalTriggerGeometryBase::getNeighborsFromTriggerCell().

Referenced by clusterizeDRNN(), and removeUnconnectedTCinCluster().

370  {
371 
372  const auto neighbors = triggerGeometry.getNeighborsFromTriggerCell( detIDa );
373 
374  if( neighbors.find( detIDb ) != neighbors.end() ) return true;
375 
376  return false;
377 
378 }
virtual geom_set getNeighborsFromTriggerCell(const unsigned trigger_cell_det_id) const =0
void HGCalClusteringImpl::calibratePt ( l1t::HGCalCluster cluster)
private

Definition at line 432 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().

432  {
433 
434  double calibPt=0.;
435 
436  if(applyLayerWeights_){
437 
438  unsigned layerN = triggerTools_.layerWithOffset(cluster.detId());
439 
440  if(layerWeights_.at(layerN)==0.){
441  throw cms::Exception("BadConfiguration")
442  <<"2D cluster energy forced to 0 by calibration coefficients.\n"
443  <<"The configuration should be changed. "
444  <<"Discarded layers should be defined in hgcalTriggerGeometryESProducer.TriggerGeometry.DisconnectedLayers and not with calibration coefficients = 0\n";
445  }
446 
447  calibPt = layerWeights_.at(layerN) * cluster.mipPt();
448 
449  }
450 
451  else{
452 
453  calibPt = cluster.pt() * calibSF_;
454 
455  }
456 
457  math::PtEtaPhiMLorentzVector calibP4( calibPt,
458  cluster.eta(),
459  cluster.phi(),
460  0. );
461 
462  cluster.setP4( calibP4 );
463 
464 }
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, createfilelist::int, 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 
78  double minDist = dr_;
79  int targetClu = -1;
80 
81  for(unsigned iclu=0; iclu<clustersTmp.size(); iclu++){
82 
83  if(!this->isPertinent(**tc, clustersTmp.at(iclu), dr_)) continue;
84 
85  double d = clustersTmp.at(iclu).distance(**tc);
86  if(d<minDist){
87  minDist = d;
88  targetClu = int(iclu);
89  }
90  }
91 
92  if(targetClu<0 && isSeed[itc]) clustersTmp.emplace_back( *tc );
93  else if(targetClu>=0) clustersTmp.at( targetClu ).addConstituent( *tc );
94 
95  }
96 
97  /* store clusters in the persistent collection */
98  clusters.resize(0, clustersTmp.size());
99  for( unsigned i(0); i<clustersTmp.size(); ++i ){
100  calibratePt(clustersTmp.at(i));
101  clusters.set( 0, i, clustersTmp.at(i) );
102  }
103 
104 
105 
106 }
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 272 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().

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

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

Referenced by eventSetup().

251  {
252 
253  std::array< std::vector< std::vector<edm::Ptr<l1t::HGCalTriggerCell>>>,kNSides_> reshuffledTriggerCells;
255  for(unsigned side=0; side<kNSides_; side++)
256  {
257  reshuffledTriggerCells[side].resize(layers);
258  }
259  triggerCellReshuffling( triggerCellsPtrs, reshuffledTriggerCells );
260 
261  for(unsigned iec=0; iec<kNSides_; ++iec){
262  for(unsigned il=0; il<layers; ++il){
263  NNKernel( reshuffledTriggerCells[iec][il], clusters, triggerGeometry );
264  }
265  }
266 
267 }
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 130 of file HGCalClusteringImpl.cc.

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

Referenced by eventSetup(), and NNKernel().

132 {
133 
134  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& pertinentTC = secondary_cluster.constituents();
135 
136  for(const auto& id_tc : pertinentTC){
137  main_cluster.addConstituent(id_tc.second);
138  }
139 
140 }
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 143 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().

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

381  {
382 
383  /* get the constituents and the centre of the seed tc (considered as the first of the constituents) */
384  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& constituents = cluster.constituents();
385  Basic3DVector<float> seedCentre( constituents.at(cluster.detId())->position() );
386 
387  /* distances from the seed */
388  vector<pair<edm::Ptr<l1t::HGCalTriggerCell>,float>> distances;
389  for( const auto & id_tc : constituents )
390  {
391  Basic3DVector<float> tcCentre( id_tc.second->position() );
392  float distance = ( seedCentre - tcCentre ).mag();
393  distances.push_back( pair<edm::Ptr<l1t::HGCalTriggerCell>,float>( id_tc.second, distance ) );
394  }
395 
396  /* sorting (needed in order to be sure that we are skipping any tc) */
397  /* FIXME: better sorting needed!!! */
398  std::sort( distances.begin(), distances.end(), distanceSorter );
399 
400 
401  /* checking if the tc is connected to the seed */
402  bool toRemove[constituents.size()];
403  toRemove[0] = false; // this is the seed
404  for( unsigned itc=1; itc<distances.size(); itc++ ){
405 
406  /* get the tc under study */
407  toRemove[itc] = true;
408  const edm::Ptr<l1t::HGCalTriggerCell>& tcToStudy = distances[itc].first;
409 
410  /* compare with the tc in the cluster */
411  for( unsigned itc_ref=0; itc_ref<itc; itc_ref++ ){
412  if( !toRemove[itc_ref] ) {
413  if( areTCneighbour( tcToStudy->detId(), distances.at( itc_ref ).first->detId(), triggerGeometry ) ) {
414  toRemove[itc] = false;
415  break;
416  }
417  }
418  }
419 
420  }
421 
422 
423  /* remove the unconnected TCs */
424  for( unsigned i=0; i<distances.size(); i++){
425  if( toRemove[i] ) cluster.removeConstituent( distances.at( i ).first );
426  }
427 
428 }
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 113 of file HGCalClusteringImpl.cc.

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

Referenced by clusterizeNN().

115  {
116 
117  for( const auto& tc : triggerCellsPtrs ){
118  int endcap = tc->zside() == -1 ? 0 : 1 ;
119  HGCalDetId tcDetId( tc->detId() );
120  unsigned layer = triggerTools_.layerWithOffset(tc->detId());
121 
122  reshuffledTriggerCells[endcap][layer-1].emplace_back(tc);
123 
124  }
125 
126 }
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