1 #include <unordered_set>
2 #include <unordered_map>
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")) {
29 double distXY)
const {
44 bool isSeed[triggerCellsPtrs.size()];
49 tc != triggerCellsPtrs.end();
52 isSeed[itc] = ((*tc)->mipPt() > seedThreshold) ?
true :
false;
56 std::vector<l1t::HGCalCluster> clustersTmp;
60 tc != triggerCellsPtrs.end();
76 for (
unsigned iclu = 0; iclu < clustersTmp.size(); iclu++) {
80 double d = clustersTmp.at(iclu).distance(**tc);
83 targetClu = int(iclu);
87 if (targetClu < 0 && isSeed[itc])
88 clustersTmp.emplace_back(*tc);
89 else if (targetClu >= 0)
90 clustersTmp.at(targetClu).addConstituent(*tc);
94 clusters.
resize(0, clustersTmp.size());
95 for (
unsigned i(0);
i < clustersTmp.size(); ++
i) {
97 clusters.
set(0,
i, clustersTmp.at(
i));
107 for (
const auto& tc : triggerCellsPtrs) {
108 DetId tcDetId(tc->detId());
112 reshuffledTriggerCells[
endcap][layer - 1].emplace_back(tc);
119 const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& pertinentTC = secondary_cluster.
constituents();
121 for (
const auto& id_tc : pertinentTC) {
130 std::vector<l1t::HGCalCluster> clustersTmp;
133 std::unordered_map<uint32_t, unsigned> cluNNmap;
136 for (
const auto& tc_ptr : reshuffledTriggerCells) {
146 bool createNewC2d(
true);
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);
155 cluNNmap.emplace(tc_ptr->detId(), tc_cluster_itr->second);
158 <<
"Trying to access a non-existing cluster. But it should exist...\n";
164 clustersTmp.emplace_back(tc_ptr);
165 clustersTmp.back().setValid(
true);
167 cluNNmap.emplace(tc_ptr->detId(), clustersTmp.size() - 1);
173 for (
auto& cluster1 : clustersTmp) {
175 if (!cluster1.valid())
179 std::unordered_set<uint32_t> cluTcSet;
180 for (
const auto& tc_clu1 : cluster1.constituents()) {
181 cluTcSet.insert(tc_clu1.second->detId());
183 for (
const auto neighbor : neighbors) {
184 cluTcSet.insert(neighbor);
188 for (
auto& cluster2 : clustersTmp) {
190 if (cluster1.detId() == cluster2.detId())
192 if (!cluster2.valid())
196 for (
const auto& tc_clu2 : cluster2.constituents()) {
197 if (cluTcSet.find(tc_clu2.second->detId()) != cluTcSet.end()) {
199 cluTcSet.insert(tc_clu2.second->detId());
201 for (
const auto neighbor : neighbors) {
202 cluTcSet.insert(neighbor);
204 cluster2.setValid(
false);
213 for (
auto& cluster : clustersTmp) {
214 if (!cluster.valid())
216 bool saveInCollection(
false);
217 for (
const auto& id_tc : cluster.constituents()) {
220 if (id_tc.second->mipPt() > seedThreshold) {
221 saveInCollection =
true;
225 if (saveInCollection) {
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);
242 for (
unsigned iec = 0; iec <
kNSides_; ++iec) {
243 for (
unsigned il = 0; il <
layers; ++il) {
244 NNKernel(reshuffledTriggerCells[iec][il], clusters, triggerGeometry);
253 bool isSeed[triggerCellsPtrs.size()];
254 std::vector<unsigned> seedPositions;
255 seedPositions.reserve(triggerCellsPtrs.size());
260 tc != triggerCellsPtrs.end();
265 isSeed[itc] = ((*tc)->mipPt() > seedThreshold) ?
true :
false;
267 seedPositions.push_back(itc);
270 for (
auto pos : seedPositions) {
271 if (((*tc)->position() - triggerCellsPtrs[pos]->position()).
mag() <
dr_) {
272 if (this->
areTCneighbour((*tc)->detId(), triggerCellsPtrs[pos]->detId(), triggerGeometry)) {
274 seedPositions.pop_back();
282 std::vector<l1t::HGCalCluster> clustersTmp;
285 clustersTmp.reserve(seedPositions.size());
286 for (
auto pos : seedPositions) {
287 clustersTmp.emplace_back(triggerCellsPtrs[pos]);
293 tc != triggerCellsPtrs.end();
305 std::vector<unsigned> tcPertinentClusters;
308 for (
const auto& clu : clustersTmp) {
310 tcPertinentClusters.push_back(iclu);
315 if (tcPertinentClusters.empty()) {
317 }
else if (tcPertinentClusters.size() == 1) {
318 clustersTmp.at(tcPertinentClusters.at(0)).addConstituent(*tc);
322 for (
auto clu : tcPertinentClusters) {
323 totMipt += clustersTmp.at(clu).seedMipPt();
326 for (
auto clu : tcPertinentClusters) {
327 double seedMipt = clustersTmp.at(clu).seedMipPt();
328 clustersTmp.at(clu).addConstituent(*tc,
true, seedMipt / totMipt);
334 clusters.
resize(0, clustersTmp.size());
335 for (
unsigned i(0);
i < clustersTmp.size(); ++
i) {
338 clusters.
set(0,
i, clustersTmp.at(
i));
347 if (neighbors.find(detIDb) != neighbors.end())
356 const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& constituents = cluster.
constituents();
360 vector<pair<edm::Ptr<l1t::HGCalTriggerCell>,
float>> distances;
361 for (
const auto& id_tc : constituents) {
372 bool toRemove[constituents.size()];
374 for (
unsigned itc = 1; itc < distances.size(); itc++) {
376 toRemove[itc] =
true;
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;
391 for (
unsigned i = 0;
i < distances.size();
i++) {
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";
421 cluster.
setP4(calibP4);
void addConstituent(const edm::Ptr< C > &c, bool updateCentre=true, float fraction=1.)
double distance(const l1t::HGCalTriggerCell &tc) const
double siliconSeedThreshold_
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
double pt() const final
transverse momentum
double scintillatorSeedThreshold_
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())
void calibratePt(l1t::HGCalCluster &cluster)
void mergeClusters(l1t::HGCalCluster &main_cluster, const l1t::HGCalCluster &secondary_cluster) const
virtual geom_set getNeighborsFromTriggerCell(const unsigned trigger_cell_det_id) const =0
std::vector< double > layerWeights_
constexpr std::array< uint8_t, layerIndexSize > layer
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
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)
static constexpr unsigned kNSides_
bool isPertinent(const l1t::HGCalTriggerCell &tc, const l1t::HGCalCluster &clu, double distXY) 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)
void set(int bx, unsigned i, const T &object)
Log< level::Info, false > LogInfo
HGCalTriggerTools triggerTools_
HGCalClusteringImpl(const edm::ParameterSet &conf)
void resize(int bx, unsigned size)
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)
bool areTCneighbour(uint32_t detIDa, uint32_t detIDb, const HGCalTriggerGeometryBase &triggerGeometry)
static int position[264][3]
double scintillatorTriggerCellThreshold_
double siliconTriggerCellThreshold_
double phi() const final
momentum azimuthal angle
void setP4(const LorentzVector &p4) final
set 4-momentum
void push_back(int bx, T object)
const std::unordered_map< uint32_t, edm::Ptr< C > > & constituents() const
double eta() const final
momentum pseudorapidity