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"))
36 if( (tcDetId.layer() != cluDetId.layer()) ||
37 (tcDetId.subdetId() != cluDetId.subdetId()) ||
38 (tcDetId.zside() != cluDetId.zside()) ){
53 bool isSeed[triggerCellsPtrs.size()];
59 isSeed[itc] = ( (*tc)->mipPt() >
seedThreshold) ?
true :
false;
63 std::vector<l1t::HGCalCluster> clustersTmp;
81 for(
unsigned iclu=0; iclu<clustersTmp.size(); iclu++){
85 double d = clustersTmp.at(iclu).distance(**tc);
88 targetClu =
int(iclu);
92 if(targetClu<0 && isSeed[itc]) clustersTmp.emplace_back( *tc );
93 else if(targetClu>=0) clustersTmp.at( targetClu ).addConstituent( *tc );
98 clusters.
resize(0, clustersTmp.size());
99 for(
unsigned i(0);
i<clustersTmp.size(); ++
i ){
101 clusters.
set( 0,
i, clustersTmp.at(
i) );
117 for(
const auto& tc : triggerCellsPtrs ){
118 int endcap = tc->zside() == -1 ? 0 : 1 ;
122 reshuffledTriggerCells[
endcap][layer-1].emplace_back(tc);
134 const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& pertinentTC = secondary_cluster.
constituents();
136 for(
const auto& id_tc : pertinentTC){
149 std::vector<l1t::HGCalCluster> clustersTmp;
152 std::unordered_map<uint32_t, unsigned> cluNNmap;
155 for(
const auto& tc_ptr : reshuffledTriggerCells ){
164 bool createNewC2d(
true);
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);
173 cluNNmap.emplace(tc_ptr->detId(), tc_cluster_itr->second);
177 <<
"Trying to access a non-existing cluster. But it should exist...\n";
183 clustersTmp.emplace_back(tc_ptr);
184 clustersTmp.back().setValid(
true);
186 cluNNmap.emplace(tc_ptr->detId(), clustersTmp.size()-1);
192 for(
auto& cluster1 : clustersTmp){
194 if( !cluster1.valid() )
continue;
197 std::unordered_set<uint32_t> cluTcSet;
198 for(
const auto& tc_clu1 : cluster1.constituents()){
199 cluTcSet.insert( tc_clu1.second->detId() );
201 for(
const auto neighbor : neighbors){
202 cluTcSet.insert( neighbor );
206 for(
auto& cluster2 : clustersTmp ){
208 if( cluster1.detId()==cluster2.detId() )
continue;
209 if( !cluster2.valid() )
continue;
212 for(
const auto& tc_clu2 : cluster2.constituents()){
213 if( cluTcSet.find(tc_clu2.second->detId())!=cluTcSet.end() ){
215 cluTcSet.insert( tc_clu2.second->detId() );
217 for(
const auto neighbor : neighbors){
218 cluTcSet.insert( neighbor );
220 cluster2.setValid(
false);
229 for(
auto& cluster : clustersTmp ){
230 if( !cluster.valid() )
continue;
231 bool saveInCollection(
false);
232 for(
const auto& id_tc : cluster.constituents() ){
236 saveInCollection =
true;
240 if(saveInCollection){
253 std::array< std::vector< std::vector<edm::Ptr<l1t::HGCalTriggerCell>>>,
kNSides_> reshuffledTriggerCells;
255 for(
unsigned side=0; side<
kNSides_; side++)
257 reshuffledTriggerCells[side].resize(layers);
261 for(
unsigned iec=0; iec<
kNSides_; ++iec){
262 for(
unsigned il=0; il<
layers; ++il){
263 NNKernel( reshuffledTriggerCells[iec][il], clusters, triggerGeometry );
277 bool isSeed[triggerCellsPtrs.size()];
278 std::vector<unsigned> seedPositions;
279 seedPositions.reserve( triggerCellsPtrs.size() );
288 isSeed[itc] = ( (*tc)->mipPt() >
seedThreshold) ?
true :
false;
291 seedPositions.push_back( itc );
294 for(
auto pos : seedPositions ){
295 if( ( (*tc)->position() - triggerCellsPtrs[
pos]->position() ).
mag() <
dr_ ){
296 if( this->
areTCneighbour( (*tc)->detId(), triggerCellsPtrs[
pos]->detId(), triggerGeometry ) )
299 seedPositions.pop_back();
308 std::vector<l1t::HGCalCluster> clustersTmp;
311 for(
auto pos : seedPositions ) {
312 clustersTmp.emplace_back( triggerCellsPtrs[
pos] );
323 if( (*tc)->mipPt() <
threshold )
continue;
324 if( isSeed[itc] )
continue;
327 std::vector<unsigned> tcPertinentClusters;
330 for (
const auto& clu : clustersTmp ) {
332 tcPertinentClusters.push_back( iclu );
337 if ( tcPertinentClusters.empty() ) {
340 else if( tcPertinentClusters.size() == 1 ) {
341 clustersTmp.at( tcPertinentClusters.at(0) ).addConstituent( *tc );
347 for(
auto clu : tcPertinentClusters ){
348 totMipt += clustersTmp.at( clu ).seedMipPt();
351 for(
auto clu : tcPertinentClusters ){
352 double seedMipt = clustersTmp.at( clu ).seedMipPt();
353 clustersTmp.at( clu ).addConstituent( *tc,
true, seedMipt/totMipt );
359 clusters.
resize(0, clustersTmp.size());
360 for(
unsigned i(0);
i<clustersTmp.size(); ++
i ){
363 clusters.
set( 0,
i, clustersTmp.at(
i) );
374 if( neighbors.find( detIDb ) != neighbors.end() )
return true;
384 const std::unordered_map<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& constituents = cluster.
constituents();
388 vector<pair<edm::Ptr<l1t::HGCalTriggerCell>,
float>> distances;
389 for(
const auto & id_tc : constituents )
402 bool toRemove[constituents.size()];
404 for(
unsigned itc=1; itc<distances.size(); itc++ ){
407 toRemove[itc] =
true;
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;
424 for(
unsigned i=0;
i<distances.size();
i++){
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";
462 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_
double eta() const final
momentum pseudorapidity
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
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())
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 pt() const final
transverse momentum
std::vector< double > layerWeights_
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)
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)
static const unsigned kNSides_
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