File indexing completed on 2024-04-06 12:20:42
0001 #include "L1Trigger/L1THGCal/interface/backend/HGCalStage1TruncationImpl_SA.h"
0002 #include <cmath>
0003
0004 unsigned HGCalStage1TruncationImplSA::run(const l1thgcfirmware::HGCalTriggerCellSACollection& tcs_in,
0005 const l1thgcfirmware::Stage1TruncationConfig& theConf,
0006 l1thgcfirmware::HGCalTriggerCellSACollection& tcs_out) const {
0007 unsigned sector120 = theConf.phiSector();
0008 std::unordered_map<unsigned, l1thgcfirmware::HGCalTriggerCellSACollection> tcs_per_bin;
0009
0010
0011 bool do_truncate = theConf.doTruncate();
0012 double rozmin = theConf.rozMin();
0013 double rozmax = theConf.rozMax();
0014 unsigned rozbins = theConf.rozBins();
0015 const std::vector<unsigned>& maxtcsperbin = theConf.maxTcsPerBin();
0016 const std::vector<double>& phiedges = theConf.phiEdges();
0017
0018
0019 for (const auto& tc : tcs_in) {
0020 double x = tc.x();
0021 double y = tc.y();
0022 double z = tc.z();
0023 double roverz = std::sqrt(x * x + y * y) / std::abs(z);
0024 unsigned roverzbin = rozBin(roverz, rozmin, rozmax, rozbins);
0025
0026 double phi = rotatedphi(x, y, z, sector120);
0027 int phibin = phiBin(roverzbin, phi, phiedges);
0028 if (phibin < 0)
0029 return 1;
0030 unsigned packed_bin = packBin(roverzbin, phibin);
0031
0032 tcs_per_bin[packed_bin].push_back(tc);
0033 }
0034
0035 for (auto& bin_tcs : tcs_per_bin) {
0036 std::sort(bin_tcs.second.begin(),
0037 bin_tcs.second.end(),
0038 [](const l1thgcfirmware::HGCalTriggerCell& a, const l1thgcfirmware::HGCalTriggerCell& b) -> bool {
0039 return a.mipPt() > b.mipPt();
0040 });
0041
0042 unsigned roverzbin = 0;
0043 unsigned phibin = 0;
0044 unpackBin(bin_tcs.first, roverzbin, phibin);
0045 if (roverzbin >= maxtcsperbin.size())
0046 return 1;
0047
0048 unsigned max_tc = maxtcsperbin[roverzbin];
0049 if (do_truncate && bin_tcs.second.size() > max_tc) {
0050 bin_tcs.second.resize(max_tc);
0051 }
0052
0053 for (const auto& tc : bin_tcs.second) {
0054 tcs_out.push_back(tc);
0055 }
0056 }
0057
0058 return 0;
0059 }
0060
0061 unsigned HGCalStage1TruncationImplSA::packBin(unsigned roverzbin, unsigned phibin) const {
0062 unsigned packed_bin = 0;
0063 packed_bin |= ((roverzbin & mask_roz_) << offset_roz_);
0064 packed_bin |= (phibin & mask_phi_);
0065 return packed_bin;
0066 }
0067
0068 void HGCalStage1TruncationImplSA::unpackBin(unsigned packedbin, unsigned& roverzbin, unsigned& phibin) const {
0069 roverzbin = ((packedbin >> offset_roz_) & mask_roz_);
0070 phibin = (packedbin & mask_phi_);
0071 }
0072
0073 int HGCalStage1TruncationImplSA::phiBin(unsigned roverzbin, double phi, const std::vector<double>& phiedges) const {
0074 int phi_bin = 0;
0075 if (roverzbin >= phiedges.size())
0076 return -1;
0077 double phi_edge = phiedges[roverzbin];
0078 if (phi > phi_edge)
0079 phi_bin = 1;
0080 return phi_bin;
0081 }
0082
0083 double HGCalStage1TruncationImplSA::rotatedphi(double x, double y, double z, int sector) const {
0084 if (z > 0)
0085 x = -x;
0086 double phi = std::atan2(y, x);
0087
0088 if (sector == 1) {
0089 if (phi < M_PI and phi > 0)
0090 phi = phi - (2. * M_PI / 3.);
0091 else
0092 phi = phi + (4. * M_PI / 3.);
0093 } else if (sector == 2) {
0094 phi = phi + (2. * M_PI / 3.);
0095 }
0096 return phi;
0097 }
0098
0099 unsigned HGCalStage1TruncationImplSA::rozBin(double roverz, double rozmin, double rozmax, unsigned rozbins) const {
0100 constexpr double margin = 1.001;
0101 double roz_bin_size = (rozbins > 0 ? (rozmax - rozmin) * margin / double(rozbins) : 0.);
0102 unsigned roverzbin = 0;
0103 if (roz_bin_size > 0.) {
0104 roverz -= rozmin;
0105 roverz = std::clamp(roverz, 0., rozmax - rozmin);
0106 roverzbin = unsigned(roverz / roz_bin_size);
0107 }
0108
0109 return roverzbin;
0110 }