CMS 3D CMS Logo

L1EGammaEEProducer.cc
Go to the documentation of this file.
4 
10 
11 namespace l1tp2 {
12  // we sort the clusters in pt
14  return cl1->pt() > cl2->pt();
15  }
16 }; // namespace l1tp2
17 
19  static float constexpr eta_min = 1.;
20  static float constexpr eta_max = 4.;
21  static unsigned constexpr n_eta_bins = 150;
22  int eta_bin = floor((std::abs(cl->eta()) - eta_min) / ((eta_max - eta_min) / n_eta_bins));
23  if (cl->eta() < 0)
24  return -1 * eta_bin; // bin 0 doesn't exist
25  return eta_bin;
26 }
27 
29  static constexpr float phi_min = -M_PI;
30  static constexpr float phi_max = M_PI;
31  static constexpr unsigned n_phi_bins = 63;
32  return floor(std::abs(reco::deltaPhi(cl->phi(), phi_min)) / ((phi_max - phi_min) / n_phi_bins));
33 }
34 
35 pair<int, int> get_eta_phi_bin(const l1t::HGCalMulticluster *cl) { return std::make_pair(etaBin(cl), get_phi_bin(cl)); }
36 
38 public:
39  explicit L1EGammaEEProducer(const edm::ParameterSet &);
40 
41 private:
42  void produce(edm::Event &, const edm::EventSetup &) override;
43 
46 };
47 
49  : multiclusters_token_(
50  consumes<l1t::HGCalMulticlusterBxCollection>(iConfig.getParameter<edm::InputTag>("Multiclusters"))),
51  calibrator_(iConfig.getParameter<edm::ParameterSet>("calibrationConfig")) {
52  produces<BXVector<l1t::EGamma>>("L1EGammaCollectionBXVWithCuts");
53 }
54 
56  float minEt_ = 0;
57 
58  std::unique_ptr<BXVector<l1t::EGamma>> l1EgammaBxCollection(new l1t::EGammaBxCollection);
59 
60  // retrieve clusters 3D
62  iEvent.getByToken(multiclusters_token_, multiclusters_h);
63  const l1t::HGCalMulticlusterBxCollection &multiclusters = *multiclusters_h;
64 
65  std::vector<const l1t::HGCalMulticluster *> selected_multiclusters;
66  std::map<std::pair<int, int>, std::vector<const l1t::HGCalMulticluster *>> etaphi_bins;
67 
68  // here we loop on the TPGs
69  for (auto cl3d = multiclusters.begin(0); cl3d != multiclusters.end(0); cl3d++) {
70  if (cl3d->hwQual()) {
71  if (cl3d->et() > minEt_) {
72  int hw_quality = 1; // baseline EG ID passed
73  if (std::abs(cl3d->eta()) >= 1.52) {
74  hw_quality = 2; // baseline EG ID passed + cleanup of transition region
75  }
76 
77  float calib_factor = calibrator_.calibrationFactor(cl3d->pt(), cl3d->eta());
78  l1t::EGamma eg =
79  l1t::EGamma(reco::Candidate::PolarLorentzVector(cl3d->pt() / calib_factor, cl3d->eta(), cl3d->phi(), 0.));
80  eg.setHwQual(hw_quality);
81  eg.setHwIso(1);
82  eg.setIsoEt(-1); // just temporarily as a dummy value
83  l1EgammaBxCollection->push_back(0, eg);
84  if (hw_quality == 2) {
85  // we build the EM interpreted EG object
87  cl3d->iPt(l1t::HGCalMulticluster::EnergyInterpretation::EM), cl3d->eta(), cl3d->phi(), 0.));
88  eg_emint.setHwQual(4);
89  eg_emint.setHwIso(1);
90  eg_emint.setIsoEt(-1); // just temporarily as a dummy value
91  l1EgammaBxCollection->push_back(0, eg_emint);
92  // we also prepare for the brem recovery procedure
93  selected_multiclusters.push_back(&(*cl3d));
94  auto eta_phi_bin = get_eta_phi_bin(&(*cl3d));
95  auto bucket = etaphi_bins.find(eta_phi_bin);
96  if (bucket == etaphi_bins.end()) {
97  std::vector<const l1t::HGCalMulticluster *> vec;
98  vec.push_back(&(*cl3d));
99  etaphi_bins[eta_phi_bin] = vec;
100  } else {
101  bucket->second.push_back(&(*cl3d));
102  }
103  }
104  }
105  }
106  }
107 
108  std::sort(selected_multiclusters.begin(), selected_multiclusters.end(), l1tp2::compare_cluster_pt);
109  std::set<const l1t::HGCalMulticluster *> used_clusters;
110  for (const auto &cl3d : selected_multiclusters) {
111  if (used_clusters.find(cl3d) == used_clusters.end()) {
112  float pt = cl3d->pt();
113  // we drop the Had component of the energy
114  if (cl3d->hOverE() != -1)
115  pt = cl3d->pt() / (1 + cl3d->hOverE());
116  reco::Candidate::PolarLorentzVector mom(pt, cl3d->eta(), cl3d->phi(), 0.);
118  cl3d->iPt(l1t::HGCalMulticluster::EnergyInterpretation::EM), cl3d->eta(), cl3d->phi(), 0.);
119 
120  // this is not yet used
121  used_clusters.insert(cl3d);
122  auto eta_phi_bin = get_eta_phi_bin(cl3d);
123 
124  for (int eta_bin : {eta_phi_bin.first - 1, eta_phi_bin.first, eta_phi_bin.first + 1}) {
125  for (int phi_bin : {eta_phi_bin.second - 1, eta_phi_bin.second, eta_phi_bin.second + 1}) {
126  auto bucket = etaphi_bins.find(std::make_pair(eta_bin, phi_bin));
127  if (bucket != etaphi_bins.end()) {
128  // this bucket is not empty
129  for (const auto &other_cl_ptr : bucket->second) {
130  if (used_clusters.find(other_cl_ptr) == used_clusters.end()) {
131  if (std::abs(other_cl_ptr->eta() - cl3d->eta()) < 0.02) {
132  if (std::abs(reco::deltaPhi(other_cl_ptr->phi(), cl3d->phi())) < 0.1) {
133  float pt_other = other_cl_ptr->pt();
134  if (other_cl_ptr->hOverE() != -1)
135  pt_other = other_cl_ptr->pt() / (1 + other_cl_ptr->hOverE());
136  mom += reco::Candidate::PolarLorentzVector(pt_other, other_cl_ptr->eta(), other_cl_ptr->phi(), 0.);
138  other_cl_ptr->iPt(l1t::HGCalMulticluster::EnergyInterpretation::EM),
139  other_cl_ptr->eta(),
140  other_cl_ptr->phi(),
141  0.);
142  used_clusters.insert(other_cl_ptr);
143  }
144  }
145  }
146  }
147  }
148  }
149  }
150  float calib_factor = calibrator_.calibrationFactor(mom.pt(), mom.eta());
151  l1t::EGamma eg =
152  l1t::EGamma(reco::Candidate::PolarLorentzVector(mom.pt() / calib_factor, mom.eta(), mom.phi(), 0.));
153  eg.setHwQual(3);
154  eg.setHwIso(1);
155  l1EgammaBxCollection->push_back(0, eg);
156 
157  l1t::EGamma eg_emint_brec =
158  l1t::EGamma(reco::Candidate::PolarLorentzVector(mom_eint.pt(), mom_eint.eta(), mom_eint.phi(), 0.));
159  eg_emint_brec.setHwQual(5);
160  eg_emint_brec.setHwIso(1);
161  l1EgammaBxCollection->push_back(0, eg_emint_brec);
162  }
163  }
164 
165  iEvent.put(std::move(l1EgammaBxCollection), "L1EGammaCollectionBXVWithCuts");
166 }
167 
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
double pt() const final
transverse momentum
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void produce(edm::Event &, const edm::EventSetup &) override
void setHwQual(int qual)
Definition: L1Candidate.h:31
delete x;
Definition: CaloConfig.h:22
static constexpr int n_eta_bins
edm::EDGetToken multiclusters_token_
const_iterator begin(int bx) const
int iEvent
Definition: GenABIO.cc:224
L1EGammaEECalibrator calibrator_
pair< int, int > get_eta_phi_bin(const l1t::HGCalMulticluster *cl)
int get_phi_bin(const l1t::HGCalMulticluster *cl)
L1EGammaEEProducer(const edm::ParameterSet &)
BXVector< HGCalMulticluster > HGCalMulticlusterBxCollection
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool compare_cluster_pt(const l1t::HGCalMulticluster *cl1, const l1t::HGCalMulticluster *cl2)
#define M_PI
void setHwIso(int iso)
Definition: L1Candidate.h:32
void setIsoEt(short int iso)
Definition: EGamma.cc:34
const_iterator end(int bx) const
HLT enums.
float calibrationFactor(const float &pt, const float &eta) const
int etaBin(const l1t::HGCalMulticluster *cl)
def move(src, dest)
Definition: eostools.py:511
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:38