Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:20:42

0001 #include "L1Trigger/L1THGCal/interface/backend/HGCalStage2ClusterDistribution.h"
0002 #include "DataFormats/Common/interface/PtrVector.h"
0003 
0004 //class constructor
0005 HGCalStage2ClusterDistribution::HGCalStage2ClusterDistribution(const edm::ParameterSet& conf)
0006     : roz_min_(conf.getParameter<double>("rozMin")),
0007       roz_max_(conf.getParameter<double>("rozMax")),
0008       roz_bins_(conf.getParameter<unsigned>("rozBins")),
0009       phi_edges_(conf.getParameter<std::vector<double>>("phiSectorEdges")) {
0010   if (phi_edges_.size() != roz_bins_)
0011     throw cms::Exception("HGCalStage2ClusterDistribution::BadConfig")
0012         << "Inconsistent sizes of phiSectorEdges and rozBins";
0013 
0014   constexpr double margin = 1.001;
0015   roz_bin_size_ = (roz_bins_ > 0 ? (roz_max_ - roz_min_) * margin / double(roz_bins_) : 0.);
0016 }
0017 
0018 HGCalTriggerGeometryBase::geom_set HGCalStage2ClusterDistribution::getStage2FPGAs(
0019     const unsigned stage1_fpga,
0020     const HGCalTriggerGeometryBase::geom_set& stage2_fpgas,
0021     const edm::Ptr<l1t::HGCalCluster>& tc_ptr) const {
0022   HGCalTriggerBackendDetId stage1_fpga_id(stage1_fpga);
0023   int sector120 = stage1_fpga_id.sector();
0024 
0025   const GlobalPoint& position = tc_ptr->position();
0026   double x = position.x();
0027   double y = position.y();
0028   double z = position.z();
0029   double roverz = std::sqrt(x * x + y * y) / std::abs(z);
0030   roverz = (roverz < roz_min_ ? roz_min_ : roverz);
0031   roverz = (roverz > roz_max_ ? roz_max_ : roverz);
0032   unsigned roverzbin = (roz_bin_size_ > 0. ? unsigned((roverz - roz_min_) / roz_bin_size_) : 0);
0033   double phi = rotatedphi(x, y, z, sector120);
0034   unsigned phibin = phiBin(roverzbin, phi);
0035 
0036   HGCalTriggerGeometryBase::geom_set output_fpgas;
0037 
0038   for (const auto& fpga : stage2_fpgas) {
0039     if (phibin == 0 || sector120 == HGCalTriggerBackendDetId(fpga).sector()) {
0040       output_fpgas.emplace(fpga);
0041     }
0042   }
0043 
0044   return output_fpgas;
0045 }
0046 
0047 unsigned HGCalStage2ClusterDistribution::phiBin(unsigned roverzbin, double phi) const {
0048   unsigned phi_bin = 0;
0049   if (roverzbin >= phi_edges_.size())
0050     throw cms::Exception("HGCalStage1TruncationImpl::OutOfRange") << "roverzbin index " << roverzbin << "out of range";
0051   double phi_edge = phi_edges_[roverzbin];
0052   if (phi > phi_edge)
0053     phi_bin = 1;
0054   return phi_bin;
0055 }
0056 
0057 double HGCalStage2ClusterDistribution::rotatedphi(double x, double y, double z, int sector) const {
0058   if (z > 0)
0059     x = -x;
0060   double phi = std::atan2(y, x);
0061 
0062   if (sector == 1) {
0063     if (phi < M_PI and phi > 0)
0064       phi = phi - (2. * M_PI / 3.);
0065     else
0066       phi = phi + (4. * M_PI / 3.);
0067   } else if (sector == 2) {
0068     phi = phi + (2. * M_PI / 3.);
0069   }
0070   return phi;
0071 }