Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Stage 1 sorting and truncation per r/z bin and FPGA #484

25 changes: 25 additions & 0 deletions L1Trigger/L1THGCal/interface/backend/HGCalBackendStage1Processor.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,25 @@
#ifndef __L1Trigger_L1THGCal_HGCalBackendStage1Processor_h__
#define __L1Trigger_L1THGCal_HGCalBackendStage1Processor_h__

#include "L1Trigger/L1THGCal/interface/HGCalProcessorBase.h"

#include "DataFormats/L1THGCal/interface/HGCalTriggerCell.h"
#include "DataFormats/L1THGCal/interface/HGCalCluster.h"

#include "L1Trigger/L1THGCal/interface/backend/HGCalStage1TruncationImpl.h"
#include "L1Trigger/L1THGCal/interface/backend/HGCalClusteringDummyImpl.h"

class HGCalBackendStage1Processor : public HGCalBackendLayer1ProcessorBase {
public:
HGCalBackendStage1Processor(const edm::ParameterSet& conf);

void run(const edm::Handle<l1t::HGCalTriggerCellBxCollection>& collHandle,
l1t::HGCalClusterBxCollection& collCluster2D,
const edm::EventSetup& es) override;

private:
std::unique_ptr<HGCalClusteringDummyImpl> clusteringDummy_;
std::unique_ptr<HGCalStage1TruncationImpl> truncation_;
};

#endif
40 changes: 40 additions & 0 deletions L1Trigger/L1THGCal/interface/backend/HGCalStage1TruncationImpl.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,40 @@
#ifndef __L1Trigger_L1THGCal_HGCalStage1TruncationImpl_h__
#define __L1Trigger_L1THGCal_HGCalStage1TruncationImpl_h__

#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "DataFormats/L1THGCal/interface/HGCalTriggerCell.h"
#include "L1Trigger/L1THGCal/interface/HGCalTriggerTools.h"

#include <vector>

class HGCalStage1TruncationImpl {
public:
HGCalStage1TruncationImpl(const edm::ParameterSet& conf);

void eventSetup(const edm::EventSetup& es) { triggerTools_.eventSetup(es); }

void run(uint32_t fpga_id,
const std::vector<edm::Ptr<l1t::HGCalTriggerCell>>& tcs_in,
std::vector<edm::Ptr<l1t::HGCalTriggerCell>>& tcs_out);

private:
HGCalTriggerTools triggerTools_;

static constexpr unsigned offset_roz_ = 1;
static constexpr unsigned mask_roz_ = 0x3f; // 6 bits, max 64 bins
static constexpr unsigned mask_phi_ = 1;

double roz_min_ = 0.;
double roz_max_ = 0.;
unsigned roz_bins_ = 42;
std::vector<unsigned> max_tcs_per_bin_;
std::vector<double> phi_edges_;
double roz_bin_size_ = 0.;

uint32_t packBin(unsigned roverzbin, unsigned phibin) const;
void unpackBin(unsigned packedbin, unsigned& roverzbin, unsigned& phibin) const;
unsigned phiBin(unsigned roverzbin, double phi) const;
double rotatedphi(double x, double y, double z, int sector) const;
};

#endif
33 changes: 33 additions & 0 deletions L1Trigger/L1THGCal/plugins/backend/HGCalBackendStage1Processor.cc
Original file line number Diff line number Diff line change
@@ -0,0 +1,33 @@
#include "L1Trigger/L1THGCal/interface/backend/HGCalBackendStage1Processor.h"

DEFINE_EDM_PLUGIN(HGCalBackendLayer1Factory, HGCalBackendStage1Processor, "HGCalBackendStage1Processor");

HGCalBackendStage1Processor::HGCalBackendStage1Processor(const edm::ParameterSet& conf)
: HGCalBackendLayer1ProcessorBase(conf) {
clusteringDummy_ = std::make_unique<HGCalClusteringDummyImpl>(conf.getParameterSet("C2d_parameters"));
truncation_ = std::make_unique<HGCalStage1TruncationImpl>(conf.getParameterSet("truncation_parameters"));
}

void HGCalBackendStage1Processor::run(const edm::Handle<l1t::HGCalTriggerCellBxCollection>& collHandle,
l1t::HGCalClusterBxCollection& collCluster2D,
const edm::EventSetup& es) {
if (clusteringDummy_)
clusteringDummy_->eventSetup(es);
if (truncation_)
truncation_->eventSetup(es);

std::unordered_map<uint32_t, std::vector<edm::Ptr<l1t::HGCalTriggerCell>>> tcs_per_fpga;

for (unsigned i = 0; i < collHandle->size(); ++i) {
edm::Ptr<l1t::HGCalTriggerCell> tc_ptr(collHandle, i);
uint32_t module = geometry_->getModuleFromTriggerCell(tc_ptr->detId());
uint32_t fpga = geometry_->getStage1FpgaFromModule(module);
tcs_per_fpga[fpga].push_back(tc_ptr);
}

std::vector<edm::Ptr<l1t::HGCalTriggerCell>> truncated_tcs;
for (auto& fpga_tcs : tcs_per_fpga) {
truncation_->run(fpga_tcs.first, fpga_tcs.second, truncated_tcs);
}
clusteringDummy_->clusterizeDummy(truncated_tcs, collCluster2D);
}
9 changes: 9 additions & 0 deletions L1Trigger/L1THGCal/python/customNewProcessors.py
Original file line number Diff line number Diff line change
@@ -1,4 +1,13 @@
import FWCore.ParameterSet.Config as cms
from L1Trigger.L1THGCal.hgcalBackEndLayer1Producer_cfi import stage1truncation_proc


def custom_stage1_truncation(process):
parameters = stage1truncation_proc.clone()
process.hgcalBackEndLayer1Producer.ProcessorParameters = parameters
process.hgcalBackEndLayer2Producer.InputCluster = cms.InputTag('hgcalBackEndLayer1Producer:HGCalBackendStage1Processor')
process.hgcalTowerProducer.InputTriggerCells = cms.InputTag('hgcalBackEndLayer1Producer:HGCalBackendStage1Processor')
return process

def custom_clustering_standalone(process):
process.hgcalBackEndLayer2Producer.ProcessorParameters.ProcessorName = cms.string('HGCalBackendLayer2Processor3DClusteringSA')
Expand Down
23 changes: 23 additions & 0 deletions L1Trigger/L1THGCal/python/hgcalBackEndLayer1Producer_cfi.py
Original file line number Diff line number Diff line change
Expand Up @@ -57,6 +57,29 @@
dR_cluster=cms.double(6.),
)

ntcs_72links = [ 1, 4, 13, 13, 10, 10, 8, 8, 8, 7, 7, 6, 6, 6, 6, 6, 5, 5, 5, 5, 5, 5, 5, 4, 4, 4, 4, 4, 4, 4, 4, 4, 3, 2, 2, 2, 2, 2, 2, 2, 2, 1]
ntcs_120links = [ 2, 7, 27, 24, 19, 17, 16, 15, 14, 14, 13, 13, 13, 12, 12, 12, 11, 11, 11, 10, 10, 10, 10, 10, 9, 9, 10, 9, 9, 9, 8, 8, 7, 5, 3, 3, 3, 3, 3, 3, 3, 3]

phi_edges = [0.98901991,0.72722052,0.6981317,0.87266463,0.93084227,0.90175345,
0.87266463,0.90175345,0.95993109,0.95993109,0.93084227,0.93084227,
0.95993109,0.98901991,0.95993109,0.95993109,0.95993109,0.98901991,
0.98901991,0.95993109,0.95993109,0.98901991,0.98901991,0.98901991,
0.98901991,0.98901991,1.01810873,0.98901991,0.98901991,0.98901991,
0.98901991,0.98901991,0.98901991,0.98901991,1.04719755,1.04719755,
1.04719755,1.04719755,1.01810873,1.04719755,1.01810873,1.01810873]

truncation_params = cms.PSet(rozMin=cms.double(0.07587128),
rozMax=cms.double(0.55508006),
rozBins=cms.uint32(42),
maxTcsPerBin=cms.vuint32(ntcs_72links),
phiSectorEdges=cms.vdouble(phi_edges)
)


stage1truncation_proc = cms.PSet(ProcessorName = cms.string('HGCalBackendStage1Processor'),
C2d_parameters = dummy_C2d_params.clone(),
truncation_parameters = truncation_params.clone()
)

be_proc = cms.PSet(ProcessorName = cms.string('HGCalBackendLayer1Processor2DClustering'),
C2d_parameters = dummy_C2d_params.clone()
Expand Down
102 changes: 102 additions & 0 deletions L1Trigger/L1THGCal/src/backend/HGCalStage1TruncationImpl.cc
Original file line number Diff line number Diff line change
@@ -0,0 +1,102 @@
#include "L1Trigger/L1THGCal/interface/backend/HGCalStage1TruncationImpl.h"
#include "L1Trigger/L1THGCal/interface/HGCalBackendDetId.h"
#include <cmath>

HGCalStage1TruncationImpl::HGCalStage1TruncationImpl(const edm::ParameterSet& conf)
: roz_min_(conf.getParameter<double>("rozMin")),
roz_max_(conf.getParameter<double>("rozMax")),
roz_bins_(conf.getParameter<unsigned>("rozBins")),
max_tcs_per_bin_(conf.getParameter<std::vector<unsigned>>("maxTcsPerBin")),
phi_edges_(conf.getParameter<std::vector<double>>("phiSectorEdges")) {
if (max_tcs_per_bin_.size() != roz_bins_)
throw cms::Exception("HGCalStage1TruncationImpl::BadConfig") << "Inconsistent sizes of maxTcsPerBin and rozBins";
if (phi_edges_.size() != roz_bins_)
throw cms::Exception("HGCalStage1TruncationImpl::BadConfig") << "Inconsistent sizes of phiSectorEdges and rozBins";

constexpr double margin = 1.001;
roz_bin_size_ = (roz_bins_ > 0 ? (roz_max_ - roz_min_) * margin / double(roz_bins_) : 0.);
}

void HGCalStage1TruncationImpl::run(uint32_t fpga_id,
const std::vector<edm::Ptr<l1t::HGCalTriggerCell>>& tcs_in,
std::vector<edm::Ptr<l1t::HGCalTriggerCell>>& tcs_out) {
unsigned sector120 = HGCalBackendDetId(fpga_id).sector();
std::unordered_map<unsigned, std::vector<edm::Ptr<l1t::HGCalTriggerCell>>> tcs_per_bin;

// group TCs per (r/z, phi) bins
for (const auto& tc : tcs_in) {
const GlobalPoint& position = tc->position();
double x = position.x();
double y = position.y();
double z = position.z();
double roverz = std::sqrt(x * x + y * y) / std::abs(z);
roverz = (roverz < roz_min_ ? roz_min_ : roverz);
roverz = (roverz > roz_max_ ? roz_max_ : roverz);
unsigned roverzbin = (roz_bin_size_ > 0. ? unsigned((roverz - roz_min_) / roz_bin_size_) : 0);
double phi = rotatedphi(x, y, z, sector120);
unsigned phibin = phiBin(roverzbin, phi);
unsigned packed_bin = packBin(roverzbin, phibin);

tcs_per_bin[packed_bin].push_back(tc);
}
// apply sorting and trunction in each (r/z, phi) bin
for (auto& bin_tcs : tcs_per_bin) {
std::sort(bin_tcs.second.begin(),
bin_tcs.second.end(),
[](const edm::Ptr<l1t::HGCalTriggerCell>& a, const edm::Ptr<l1t::HGCalTriggerCell>& b) -> bool {
return a->mipPt() > b->mipPt();
});

unsigned roverzbin = 0;
unsigned phibin = 0;
unpackBin(bin_tcs.first, roverzbin, phibin);
if (roverzbin >= max_tcs_per_bin_.size())
throw cms::Exception("HGCalStage1TruncationImpl::OutOfRange")
<< "roverzbin index " << roverzbin << "out of range";
unsigned max_tc = max_tcs_per_bin_[roverzbin];
if (bin_tcs.second.size() > max_tc) {
bin_tcs.second.resize(max_tc);
}
for (const auto& tc : bin_tcs.second) {
tcs_out.push_back(tc);
}
}
}

unsigned HGCalStage1TruncationImpl::packBin(unsigned roverzbin, unsigned phibin) const {
unsigned packed_bin = 0;
packed_bin |= ((roverzbin & mask_roz_) << offset_roz_);
packed_bin |= (phibin & mask_phi_);
return packed_bin;
}

void HGCalStage1TruncationImpl::unpackBin(unsigned packedbin, unsigned& roverzbin, unsigned& phibin) const {
roverzbin = ((packedbin >> offset_roz_) & mask_roz_);
phibin = (packedbin & mask_phi_);
}

unsigned HGCalStage1TruncationImpl::phiBin(unsigned roverzbin, double phi) const {
unsigned phi_bin = 0;
if (roverzbin >= phi_edges_.size())
throw cms::Exception("HGCalStage1TruncationImpl::OutOfRange") << "roverzbin index " << roverzbin << "out of range";
double phi_edge = phi_edges_[roverzbin];
if (phi > phi_edge)
phi_bin = 1;
return phi_bin;
}

double HGCalStage1TruncationImpl::rotatedphi(double x, double y, double z, int sector) const {
if (z > 0)
x = -x;
double phi = std::atan2(y, x);

if (sector == 1) {
if (phi < M_PI and phi > 0)
phi = phi - (2. * M_PI / 3.);
else
phi = phi + (4. * M_PI / 3.);
} else if (sector == 2) {
phi = phi + (2. * M_PI / 3.);
}
return phi;
}
7 changes: 7 additions & 0 deletions L1Trigger/L1THGCalUtilities/python/customNtuples.py
Original file line number Diff line number Diff line change
Expand Up @@ -8,6 +8,13 @@ def custom_ntuples_V9(process):
ntuple.bhSimHits = cms.InputTag('g4SimHits:HGCHitsHEback')
return process

def custom_ntuples_stage1_truncation(process):
ntuples = process.hgcalTriggerNtuplizer.Ntuples
for ntuple in ntuples:
if ntuple.NtupleName=='HGCalTriggerNtupleHGCClusters':
ntuple.Clusters = cms.InputTag('hgcalBackEndLayer1Producer:HGCalBackendStage1Processor')
return process

def custom_ntuples_standalone_clustering(process):
ntuples = process.hgcalTriggerNtuplizer.Ntuples
for ntuple in ntuples:
Expand Down
Original file line number Diff line number Diff line change
@@ -0,0 +1,96 @@
import FWCore.ParameterSet.Config as cms

from Configuration.Eras.Era_Phase2C9_cff import Phase2C9
process = cms.Process('DIGI',Phase2C9)

# import of standard configurations
process.load('Configuration.StandardSequences.Services_cff')
process.load('SimGeneral.HepPDTESSource.pythiapdt_cfi')
process.load('FWCore.MessageService.MessageLogger_cfi')
process.load('Configuration.EventContent.EventContent_cff')
process.load('SimGeneral.MixingModule.mixNoPU_cfi')
process.load('Configuration.Geometry.GeometryExtended2026D49Reco_cff')
process.load('Configuration.Geometry.GeometryExtended2026D49_cff')
process.load('Configuration.StandardSequences.MagneticField_cff')
process.load('Configuration.StandardSequences.Generator_cff')
process.load('IOMC.EventVertexGenerators.VtxSmearedHLLHC14TeV_cfi')
process.load('GeneratorInterface.Core.genFilterSummary_cff')
process.load('Configuration.StandardSequences.SimIdeal_cff')
process.load('Configuration.StandardSequences.Digi_cff')
process.load('Configuration.StandardSequences.SimL1Emulator_cff')
process.load('Configuration.StandardSequences.DigiToRaw_cff')
process.load('Configuration.StandardSequences.EndOfProcess_cff')
process.load('Configuration.StandardSequences.FrontierConditions_GlobalTag_cff')


process.maxEvents = cms.untracked.PSet(
input = cms.untracked.int32(5)
)

# Input source
process.source = cms.Source("PoolSource",
fileNames = cms.untracked.vstring('/store/mc/Phase2HLTTDRWinter20DIGI/SingleElectron_PT2to200/GEN-SIM-DIGI-RAW/PU200_110X_mcRun4_realistic_v3_ext2-v2/40000/00582F93-5A2A-5847-8162-D81EE503500F.root'),
inputCommands=cms.untracked.vstring(
'keep *',
'drop l1tEMTFHit2016Extras_simEmtfDigis_CSC_HLT',
'drop l1tEMTFHit2016Extras_simEmtfDigis_RPC_HLT',
'drop l1tEMTFHit2016s_simEmtfDigis__HLT',
'drop l1tEMTFTrack2016Extras_simEmtfDigis__HLT',
'drop l1tEMTFTrack2016s_simEmtfDigis__HLT',
'drop FTLClusteredmNewDetSetVector_mtdClusters_FTLBarrel_RECO',
'drop FTLClusteredmNewDetSetVector_mtdClusters_FTLEndcap_RECO',
'drop MTDTrackingRecHitedmNewDetSetVector_mtdTrackingRecHits__RECO',
'drop BTLDetIdBTLSampleFTLDataFrameTsSorted_mix_FTLBarrel_HLT',
'drop ETLDetIdETLSampleFTLDataFrameTsSorted_mix_FTLEndcap_HLT',
)
)

process.options = cms.untracked.PSet(

)

# Production Info
process.configurationMetadata = cms.untracked.PSet(
version = cms.untracked.string('$Revision: 1.20 $'),
annotation = cms.untracked.string('SingleElectronPt10_cfi nevts:10'),
name = cms.untracked.string('Applications')
)

# Output definition
process.TFileService = cms.Service(
"TFileService",
fileName = cms.string("ntuple.root")
)

# Other statements
from Configuration.AlCa.GlobalTag import GlobalTag
process.GlobalTag = GlobalTag(process.GlobalTag, 'auto:phase2_realistic_T15', '')

# load HGCAL TPG simulation
process.load('L1Trigger.L1THGCal.hgcalTriggerPrimitives_cff')

# Use new Stage 1 processor
from L1Trigger.L1THGCal.customNewProcessors import custom_stage1_truncation
process = custom_stage1_truncation(process)

# Switch to latest trigger geometry containing information on links mapping
from L1Trigger.L1THGCal.customTriggerGeometry import custom_geometry_decentralized_V11
process = custom_geometry_decentralized_V11(process, links='signaldriven', implementation=2)

process.hgcl1tpg_step = cms.Path(process.hgcalTriggerPrimitives)


# load ntuplizer and custom to use collections from Stag1 truncation processor
process.load('L1Trigger.L1THGCalUtilities.hgcalTriggerNtuples_cff')
from L1Trigger.L1THGCalUtilities.customNtuples import custom_ntuples_stage1_truncation
process = custom_ntuples_stage1_truncation(process)
process.ntuple_step = cms.Path(process.hgcalTriggerNtuples)

# Schedule definition
process.schedule = cms.Schedule(process.hgcl1tpg_step, process.ntuple_step)

# Add early deletion of temporary data products to reduce peak memory need
from Configuration.StandardSequences.earlyDeleteSettings_cff import customiseEarlyDelete
process = customiseEarlyDelete(process)
# End adding early deletion