CMS 3D CMS Logo

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

#include <HGCalMulticlusteringImpl.h>

Public Member Functions

void clusterizeDBSCAN (const std::vector< edm::Ptr< l1t::HGCalCluster >> &clustersPtr, l1t::HGCalMulticlusterBxCollection &multiclusters, const HGCalTriggerGeometryBase &triggerGeometry)
 
void clusterizeDR (const std::vector< edm::Ptr< l1t::HGCalCluster >> &clustersPtr, l1t::HGCalMulticlusterBxCollection &multiclusters, const HGCalTriggerGeometryBase &triggerGeometry)
 
 HGCalMulticlusteringImpl (const edm::ParameterSet &conf)
 
bool isPertinent (const l1t::HGCalCluster &clu, const l1t::HGCalMulticluster &mclu, double dR) const
 
void setGeometry (const HGCalTriggerGeometryBase *const geom)
 

Private Member Functions

void finalizeClusters (std::vector< l1t::HGCalMulticluster > &, l1t::HGCalMulticlusterBxCollection &, const HGCalTriggerGeometryBase &)
 
void findNeighbor (const std::vector< std::pair< unsigned int, double >> &rankedList, unsigned int searchInd, const std::vector< edm::Ptr< l1t::HGCalCluster >> &clustersPtr, std::vector< unsigned int > &neigbors)
 

Private Attributes

double distDbscan_ = 0.005
 
double dr_
 
std::unique_ptr< HGCalTriggerClusterIdentificationBaseid_
 
unsigned minNDbscan_ = 3
 
std::string multiclusterAlgoType_
 
double ptC3dThreshold_
 
HGCalShowerShape shape_
 
HGCalTriggerTools triggerTools_
 

Detailed Description

Definition at line 14 of file HGCalMulticlusteringImpl.h.

Constructor & Destructor Documentation

◆ HGCalMulticlusteringImpl()

HGCalMulticlusteringImpl::HGCalMulticlusteringImpl ( const edm::ParameterSet conf)

Definition at line 5 of file HGCalMulticlusteringImpl.cc.

References distDbscan_, dr_, get, edm::ParameterSet::getParameter(), id_, minNDbscan_, multiclusterAlgoType_, and ptC3dThreshold_.

6  : dr_(conf.getParameter<double>("dR_multicluster")),
7  ptC3dThreshold_(conf.getParameter<double>("minPt_multicluster")),
8  multiclusterAlgoType_(conf.getParameter<string>("type_multicluster")),
9  distDbscan_(conf.getParameter<double>("dist_dbscan_multicluster")),
10  minNDbscan_(conf.getParameter<unsigned>("minN_dbscan_multicluster")) {
11  edm::LogInfo("HGCalMulticlusterParameters") << "Multicluster dR for Near Neighbour search: " << dr_;
12  edm::LogInfo("HGCalMulticlusterParameters") << "Multicluster minimum transverse-momentum: " << ptC3dThreshold_;
13  edm::LogInfo("HGCalMulticlusterParameters") << "Multicluster DBSCAN Clustering distance: " << distDbscan_;
14  edm::LogInfo("HGCalMulticlusterParameters") << "Multicluster clustering min number of subclusters: " << minNDbscan_;
15  edm::LogInfo("HGCalMulticlusterParameters")
16  << "Multicluster type of multiclustering algortihm: " << multiclusterAlgoType_;
17  id_ = std::unique_ptr<HGCalTriggerClusterIdentificationBase>{
18  HGCalTriggerClusterIdentificationFactory::get()->create("HGCalTriggerClusterIdentificationBDT")};
19  id_->initialize(conf.getParameter<edm::ParameterSet>("EGIdentification"));
20 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
std::unique_ptr< HGCalTriggerClusterIdentificationBase > id_
Log< level::Info, false > LogInfo
#define get

Member Function Documentation

◆ clusterizeDBSCAN()

void HGCalMulticlusteringImpl::clusterizeDBSCAN ( const std::vector< edm::Ptr< l1t::HGCalCluster >> &  clustersPtr,
l1t::HGCalMulticlusterBxCollection multiclusters,
const HGCalTriggerGeometryBase triggerGeometry 
)

Definition at line 103 of file HGCalMulticlusteringImpl.cc.

References finalizeClusters(), findNeighbor(), minNDbscan_, eostools::move(), jetUpdater_cfi::sort, triggerTools_, trackerHitRTTI::vector, class-composition::visited, and HGCalTriggerTools::zside().

105  {
106  std::vector<l1t::HGCalMulticluster> multiclustersTmp;
107  l1t::HGCalMulticluster mcluTmp;
108  std::vector<bool> visited(clustersPtrs.size(), false);
109  std::vector<bool> merged(clustersPtrs.size(), false);
110  std::vector<std::pair<unsigned int, double>> rankedList;
111  rankedList.reserve(clustersPtrs.size());
112  std::vector<std::vector<unsigned int>> neighborList;
113  neighborList.reserve(clustersPtrs.size());
114 
115  int iclu = 0, imclu = 0, neighNo;
116  double dist = 0.;
117 
118  for (std::vector<edm::Ptr<l1t::HGCalCluster>>::const_iterator clu = clustersPtrs.begin(); clu != clustersPtrs.end();
119  ++clu, ++iclu) {
120  dist = (*clu)->centreProj().mag() * triggerTools_.zside((*clu)->detId());
121  rankedList.push_back(std::make_pair(iclu, dist));
122  }
123  iclu = 0;
124  std::sort(rankedList.begin(), rankedList.end(), [](auto& left, auto& right) { return left.second < right.second; });
125 
126  for (const auto& cluRanked : rankedList) {
127  std::vector<unsigned int> neighbors;
128 
129  if (!visited.at(iclu)) {
130  visited.at(iclu) = true;
131  findNeighbor(rankedList, iclu, clustersPtrs, neighbors);
132  neighborList.push_back(std::move(neighbors));
133 
134  if (neighborList.at(iclu).size() >= minNDbscan_) {
135  multiclustersTmp.emplace_back(clustersPtrs[cluRanked.first]);
136  merged.at(iclu) = true;
137  /* dynamic range loop: range-based loop syntax cannot be employed */
138  for (unsigned int neighInd = 0; neighInd < neighborList.at(iclu).size(); neighInd++) {
139  neighNo = neighborList.at(iclu).at(neighInd);
140  /* This condition also ensures merging of clusters visited by other clusters but not merged. */
141  if (!merged.at(neighNo)) {
142  merged.at(neighNo) = true;
143  multiclustersTmp.at(imclu).addConstituent(clustersPtrs[rankedList.at(neighNo).first]);
144 
145  if (!visited.at(neighNo)) {
146  visited.at(neighNo) = true;
147  std::vector<unsigned int> secNeighbors;
148  findNeighbor(rankedList, neighNo, clustersPtrs, secNeighbors);
149 
150  if (secNeighbors.size() >= minNDbscan_) {
151  neighborList.at(iclu).insert(neighborList.at(iclu).end(), secNeighbors.begin(), secNeighbors.end());
152  }
153  }
154  }
155  }
156  imclu++;
157  }
158  } else
159  neighborList.push_back(std::move(neighbors));
160  iclu++;
161  }
162  /* making the collection of multiclusters */
163  finalizeClusters(multiclustersTmp, multiclusters, triggerGeometry);
164 }
int zside(const DetId &) const
void findNeighbor(const std::vector< std::pair< unsigned int, double >> &rankedList, unsigned int searchInd, const std::vector< edm::Ptr< l1t::HGCalCluster >> &clustersPtr, std::vector< unsigned int > &neigbors)
void finalizeClusters(std::vector< l1t::HGCalMulticluster > &, l1t::HGCalMulticlusterBxCollection &, const HGCalTriggerGeometryBase &)
def move(src, dest)
Definition: eostools.py:511

◆ clusterizeDR()

void HGCalMulticlusteringImpl::clusterizeDR ( const std::vector< edm::Ptr< l1t::HGCalCluster >> &  clustersPtr,
l1t::HGCalMulticlusterBxCollection multiclusters,
const HGCalTriggerGeometryBase triggerGeometry 
)

Definition at line 72 of file HGCalMulticlusteringImpl.cc.

References ztail::d, dr_, finalizeClusters(), createfilelist::int, isPertinent(), mag(), and trackerHitRTTI::vector.

74  {
75  std::vector<l1t::HGCalMulticluster> multiclustersTmp;
76 
77  int iclu = 0;
78  for (std::vector<edm::Ptr<l1t::HGCalCluster>>::const_iterator clu = clustersPtrs.begin(); clu != clustersPtrs.end();
79  ++clu, ++iclu) {
80  double minDist = dr_;
81  int targetMulticlu = -1;
82 
83  for (unsigned imclu = 0; imclu < multiclustersTmp.size(); imclu++) {
84  if (!this->isPertinent(**clu, multiclustersTmp.at(imclu), dr_))
85  continue;
86 
87  double d = (multiclustersTmp.at(imclu).centreProj() - (*clu)->centreProj()).mag();
88  if (d < minDist) {
89  minDist = d;
90  targetMulticlu = int(imclu);
91  }
92  }
93 
94  if (targetMulticlu < 0)
95  multiclustersTmp.emplace_back(*clu);
96  else
97  multiclustersTmp.at(targetMulticlu).addConstituent(*clu);
98  }
99 
100  /* making the collection of multiclusters */
101  finalizeClusters(multiclustersTmp, multiclusters, triggerGeometry);
102 }
d
Definition: ztail.py:151
void finalizeClusters(std::vector< l1t::HGCalMulticluster > &, l1t::HGCalMulticlusterBxCollection &, const HGCalTriggerGeometryBase &)
bool isPertinent(const l1t::HGCalCluster &clu, const l1t::HGCalMulticluster &mclu, double dR) const
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())

◆ finalizeClusters()

void HGCalMulticlusteringImpl::finalizeClusters ( std::vector< l1t::HGCalMulticluster > &  multiclusters_in,
l1t::HGCalMulticlusterBxCollection multiclusters_out,
const HGCalTriggerGeometryBase triggerGeometry 
)
private

Definition at line 166 of file HGCalMulticlusteringImpl.cc.

References bsc_activity_cfg::clusters, HGCalShowerShape::fillShapes(), id_, ptC3dThreshold_, BXVector< T >::push_back(), shape_, and TtFullHadEvtBuilder_cfi::sumPt.

Referenced by clusterizeDBSCAN(), and clusterizeDR().

168  {
169  for (auto& multicluster : multiclusters_in) {
170  // compute the eta, phi from its barycenter
171  // + pT as scalar sum of pT of constituents
172  double sumPt = 0.;
173  const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalCluster>>& clusters = multicluster.constituents();
174  for (const auto& id_cluster : clusters)
175  sumPt += id_cluster.second->pt();
176 
177  math::PtEtaPhiMLorentzVector multiclusterP4(sumPt, multicluster.centre().eta(), multicluster.centre().phi(), 0.);
178  multicluster.setP4(multiclusterP4);
179 
180  if (multicluster.pt() > ptC3dThreshold_) {
181  //compute shower shapes
182  shape_.fillShapes(multicluster, triggerGeometry);
183  // fill quality flag
184  multicluster.setHwQual(id_->decision(multicluster));
185  // fill H/E
186  multicluster.saveHOverE();
187 
188  multiclusters_out.push_back(0, multicluster);
189  }
190  }
191 }
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
void fillShapes(l1t::HGCalMulticluster &, const HGCalTriggerGeometryBase &) const
std::unique_ptr< HGCalTriggerClusterIdentificationBase > id_
void push_back(int bx, T object)

◆ findNeighbor()

void HGCalMulticlusteringImpl::findNeighbor ( const std::vector< std::pair< unsigned int, double >> &  rankedList,
unsigned int  searchInd,
const std::vector< edm::Ptr< l1t::HGCalCluster >> &  clustersPtr,
std::vector< unsigned int > &  neigbors 
)
private

Definition at line 37 of file HGCalMulticlusteringImpl.cc.

References distDbscan_, Exception, and mag().

Referenced by clusterizeDBSCAN().

40  {
41  if (clustersPtrs.size() <= searchInd || clustersPtrs.size() < rankedList.size()) {
42  throw cms::Exception("IndexOutOfBound: clustersPtrs in 'findNeighbor'");
43  }
44 
45  for (unsigned int ind = searchInd + 1;
46  ind < rankedList.size() && fabs(rankedList.at(ind).second - rankedList.at(searchInd).second) < distDbscan_;
47  ind++) {
48  if (clustersPtrs.size() <= rankedList.at(ind).first) {
49  throw cms::Exception("IndexOutOfBound: clustersPtrs in 'findNeighbor'");
50 
51  } else if (((*(clustersPtrs[rankedList.at(ind).first])).centreProj() -
52  (*(clustersPtrs[rankedList.at(searchInd).first])).centreProj())
53  .mag() < distDbscan_) {
54  neighbors.push_back(ind);
55  }
56  }
57 
58  for (unsigned int ind = 0;
59  ind < searchInd && fabs(rankedList.at(searchInd).second - rankedList.at(ind).second) < distDbscan_;
60  ind++) {
61  if (clustersPtrs.size() <= rankedList.at(ind).first) {
62  throw cms::Exception("IndexOutOfBound: clustersPtrs in 'findNeighbor'");
63 
64  } else if (((*(clustersPtrs[rankedList.at(ind).first])).centreProj() -
65  (*(clustersPtrs[rankedList.at(searchInd).first])).centreProj())
66  .mag() < distDbscan_) {
67  neighbors.push_back(ind);
68  }
69  }
70 }
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())

◆ isPertinent()

bool HGCalMulticlusteringImpl::isPertinent ( const l1t::HGCalCluster clu,
const l1t::HGCalMulticluster mclu,
double  dR 
) const

Definition at line 22 of file HGCalMulticlusteringImpl.cc.

References l1t::HGCalClusterT< C >::centreProj(), l1t::HGCalClusterT< C >::detId(), HGC3DClusterGenMatchSelector_cfi::dR, mag(), triggerTools_, and HGCalTriggerTools::zside().

Referenced by clusterizeDR().

24  {
25  DetId cluDetId(clu.detId());
26  DetId firstClusterDetId(mclu.detId());
27 
28  if (triggerTools_.zside(cluDetId) != triggerTools_.zside(firstClusterDetId)) {
29  return false;
30  }
31  if ((mclu.centreProj() - clu.centreProj()).mag() < dR) {
32  return true;
33  }
34  return false;
35 }
uint32_t detId() const
Definition: HGCalClusterT.h:92
int zside(const DetId &) const
Definition: DetId.h:17
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
const GlobalPoint & centreProj() const

◆ setGeometry()

void HGCalMulticlusteringImpl::setGeometry ( const HGCalTriggerGeometryBase *const  geom)
inline

Definition at line 18 of file HGCalMulticlusteringImpl.h.

References relativeConstraints::geom, HGCalShowerShape::setGeometry(), HGCalTriggerTools::setGeometry(), shape_, and triggerTools_.

18  {
21  }
void setGeometry(const HGCalTriggerGeometryBase *const)
void setGeometry(const HGCalTriggerGeometryBase *const geom)

Member Data Documentation

◆ distDbscan_

double HGCalMulticlusteringImpl::distDbscan_ = 0.005
private

Definition at line 45 of file HGCalMulticlusteringImpl.h.

Referenced by findNeighbor(), and HGCalMulticlusteringImpl().

◆ dr_

double HGCalMulticlusteringImpl::dr_
private

Definition at line 42 of file HGCalMulticlusteringImpl.h.

Referenced by clusterizeDR(), and HGCalMulticlusteringImpl().

◆ id_

std::unique_ptr<HGCalTriggerClusterIdentificationBase> HGCalMulticlusteringImpl::id_
private

Definition at line 50 of file HGCalMulticlusteringImpl.h.

Referenced by finalizeClusters(), and HGCalMulticlusteringImpl().

◆ minNDbscan_

unsigned HGCalMulticlusteringImpl::minNDbscan_ = 3
private

Definition at line 46 of file HGCalMulticlusteringImpl.h.

Referenced by clusterizeDBSCAN(), and HGCalMulticlusteringImpl().

◆ multiclusterAlgoType_

std::string HGCalMulticlusteringImpl::multiclusterAlgoType_
private

Definition at line 44 of file HGCalMulticlusteringImpl.h.

Referenced by HGCalMulticlusteringImpl().

◆ ptC3dThreshold_

double HGCalMulticlusteringImpl::ptC3dThreshold_
private

Definition at line 43 of file HGCalMulticlusteringImpl.h.

Referenced by finalizeClusters(), and HGCalMulticlusteringImpl().

◆ shape_

HGCalShowerShape HGCalMulticlusteringImpl::shape_
private

Definition at line 48 of file HGCalMulticlusteringImpl.h.

Referenced by finalizeClusters(), and setGeometry().

◆ triggerTools_

HGCalTriggerTools HGCalMulticlusteringImpl::triggerTools_
private

Definition at line 49 of file HGCalMulticlusteringImpl.h.

Referenced by clusterizeDBSCAN(), isPertinent(), and setGeometry().