CMS 3D CMS Logo

L1CaloJetProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: L1CaloTrigger
4 // Class: L1CaloJetProducer
5 //
16 //
17 // Original Author: Tyler Ruggles
18 // Created: Tue Aug 29 2018
19 // $Id$
20 //
21 //
22 
25 
33 
34 #include <iostream>
35 
37 //#include "DataFormats/L1TCalorimeterPhase2/interface/CaloCrystalCluster.h"
42 
45 
47 
48 // For pT calibrations
49 #include "TF1.h"
50 
51 // Run2/PhaseI output formats
54 
56 public:
57  explicit L1CaloJetProducer(const edm::ParameterSet &);
58 
59 private:
60  void produce(edm::Event &, const edm::EventSetup &) override;
61  //bool cluster_passes_base_cuts(float &cluster_pt, float &cluster_eta, float &iso, float &e2x5, float &e5x5) const;
62  int tower_diPhi(int &iPhi_1, int &iPhi_2) const;
63  int tower_diEta(int &iEta_1, int &iEta_2) const;
65  float get_hcal_calibration(float &jet_pt, float &ecal_pt, float &ecal_L1EG_jet_pt, float &jet_eta) const;
66  float get_tau_pt_calibration(float &tau_pt, float &ecal_pt, float &l1EG_pt, float &n_L1EGs, float &tau_eta) const;
67  int loose_iso_tau_wp(float &tau_pt, float &tau_iso_et, float &tau_eta) const;
68 
69  double HcalTpEtMin;
70  double EcalTpEtMin;
73  double HFTpEtMin;
77 
78  // For fetching jet calibrations
79  std::vector<double> jetPtBins;
80  std::vector<double> emFractionBinsBarrel;
81  std::vector<double> absEtaBinsBarrel;
82  std::vector<double> jetCalibrationsBarrel;
83  std::vector<double> emFractionBinsHGCal;
84  std::vector<double> absEtaBinsHGCal;
85  std::vector<double> jetCalibrationsHGCal;
86  std::vector<double> emFractionBinsHF;
87  std::vector<double> absEtaBinsHF;
88  std::vector<double> jetCalibrationsHF;
89 
90  // For fetching tau calibrations
91  std::vector<double> tauPtBins;
92  std::vector<double> tauAbsEtaBinsBarrel;
93  std::vector<double> tauCalibrationsBarrel;
94  std::vector<edm::ParameterSet> tauL1egInfoBarrel;
95  std::vector<double> tauAbsEtaBinsHGCal;
96  std::vector<double> tauCalibrationsHGCal;
97  std::vector<edm::ParameterSet> tauL1egInfoHGCal;
98 
99  // For storing jet calibrations
100  std::vector<std::vector<std::vector<double>>> calibrationsBarrel;
101  std::vector<std::vector<std::vector<double>>> calibrationsHGCal;
102  std::vector<std::vector<std::vector<double>>> calibrationsHF;
103 
104  // For storing tau calibration info
105  std::map<double, std::vector<double>> tauL1egInfoMapBarrel;
106  std::map<double, std::vector<double>> tauL1egInfoMapHGCal;
107  std::vector<double> tauL1egValuesBarrel; // To preserve ordering
108  std::vector<double> tauL1egValuesHGCal; // To preserve ordering
109  std::vector<std::vector<std::vector<std::vector<double>>>> tauPtCalibrationsBarrel;
110  std::vector<std::vector<std::vector<std::vector<double>>>> tauPtCalibrationsHGCal;
111 
112  bool debug;
115 
116  // TF1s defining tau isolation thresholds
117  TF1 isoTauBarrel = TF1("isoTauBarrelFunction", "([0] + [1]*TMath::Exp(-[2]*x))");
118  TF1 isoTauHGCal = TF1("isoTauHGCalFunction", "([0] + [1]*TMath::Exp(-[2]*x))");
119 
120  class l1CaloJetObj {
121  public:
122  bool barrelSeeded = true; // default to barrel seeded
127  //reco::Candidate::PolarLorentzVector leadingL1EG;
129  float jetClusterET = 0.;
130  float hcalJetClusterET = 0.;
131  float ecalJetClusterET = 0.;
132  float seedTowerET = 0.;
133  float leadingL1EGET = 0.;
134  float l1egJetClusterET = 0.;
135 
136  // For decay mode related checks with CaloTaus
137  std::vector<std::vector<float>> associated_l1EGs_;
138 
139  int seed_iEta = -99;
140  int seed_iPhi = -99;
141 
142  float hcal_seed = 0.;
143  //float hcal_3x3 = 0.;
144  float hcal_3x5 = 0.;
145  //float hcal_5x5 = 0.;
146  //float hcal_5x7 = 0.;
147  float hcal_7x7 = 0.;
148  //float hcal_2x3 = 0.;
149  //float hcal_2x3_1 = 0.;
150  //float hcal_2x3_2 = 0.;
151  //float hcal_2x2 = 0.;
152  //float hcal_2x2_1 = 0.;
153  //float hcal_2x2_2 = 0.;
154  //float hcal_2x2_3 = 0.;
155  //float hcal_2x2_4 = 0.;
156  float hcal_nHits = 0.;
157 
158  float ecal_seed = 0.;
159  //float ecal_3x3 = 0.;
160  float ecal_3x5 = 0.;
161  //float ecal_5x5 = 0.;
162  //float ecal_5x7 = 0.;
163  float ecal_7x7 = 0.;
164  //float ecal_2x3 = 0.;
165  //float ecal_2x3_1 = 0.;
166  //float ecal_2x3_2 = 0.;
167  //float ecal_2x2 = 0.;
168  //float ecal_2x2_1 = 0.;
169  //float ecal_2x2_2 = 0.;
170  //float ecal_2x2_3 = 0.;
171  //float ecal_2x2_4 = 0.;
172  float ecal_nHits = 0.;
173 
174  float l1eg_seed = 0.;
175  //float l1eg_3x3 = 0.;
176  float l1eg_3x5 = 0.;
177  //float l1eg_5x5 = 0.;
178  //float l1eg_5x7 = 0.;
179  float l1eg_7x7 = 0.;
180  //float l1eg_2x3 = 0.;
181  //float l1eg_2x3_1 = 0.;
182  //float l1eg_2x3_2 = 0.;
183  //float l1eg_2x2 = 0.;
184  //float l1eg_2x2_1 = 0.;
185  //float l1eg_2x2_2 = 0.;
186  //float l1eg_2x2_3 = 0.;
187  //float l1eg_2x2_4 = 0.;
188  float l1eg_nHits = 0.;
190 
191  float l1eg_nL1EGs = 0.;
196 
197  float total_seed = 0.;
198  //float total_3x3 = 0.;
199  float total_3x5 = 0.;
200  //float total_5x5 = 0.;
201  //float total_5x7 = 0.;
202  float total_7x7 = 0.;
203  //float total_2x3 = 0.;
204  //float total_2x3_1 = 0.;
205  //float total_2x3_2 = 0.;
206  //float total_2x2 = 0.;
207  //float total_2x2_1 = 0.;
208  //float total_2x2_2 = 0.;
209  //float total_2x2_3 = 0.;
210  //float total_2x2_4 = 0.;
211  float total_nHits = 0.;
212 
213  void Init() {
214  SetJetClusterP4(0., 0., 0., 0.);
215  SetHcalJetClusterP4(0., 0., 0., 0.);
216  SetEcalJetClusterP4(0., 0., 0., 0.);
217  SetSeedP4(0., 0., 0., 0.);
218  //SetLeadingL1EGP4( 0., 0., 0., 0. );
219  SetL1EGJetP4(0., 0., 0., 0.);
220  }
221 
222  void SetJetClusterP4(double pt, double eta, double phi, double mass) {
223  this->jetCluster.SetPt(pt);
224  this->jetCluster.SetEta(eta);
225  this->jetCluster.SetPhi(phi);
226  this->jetCluster.SetM(mass);
227  }
228  void SetHcalJetClusterP4(double pt, double eta, double phi, double mass) {
229  this->hcalJetCluster.SetPt(pt);
230  this->hcalJetCluster.SetEta(eta);
231  this->hcalJetCluster.SetPhi(phi);
232  this->hcalJetCluster.SetM(mass);
233  }
234  void SetEcalJetClusterP4(double pt, double eta, double phi, double mass) {
235  this->ecalJetCluster.SetPt(pt);
236  this->ecalJetCluster.SetEta(eta);
237  this->ecalJetCluster.SetPhi(phi);
238  this->ecalJetCluster.SetM(mass);
239  }
240  void SetSeedP4(double pt, double eta, double phi, double mass) {
241  this->seedTower.SetPt(pt);
242  this->seedTower.SetEta(eta);
243  this->seedTower.SetPhi(phi);
244  this->seedTower.SetM(mass);
245  }
246  //void SetLeadingL1EGP4( double pt, double eta, double phi, double mass )
247  //{
248  // this->leadingL1EG.SetPt( pt );
249  // this->leadingL1EG.SetEta( eta );
250  // this->leadingL1EG.SetPhi( phi );
251  // this->leadingL1EG.SetM( mass );
252  //}
253  void SetL1EGJetP4(double pt, double eta, double phi, double mass) {
254  this->l1egJetCluster.SetPt(pt);
255  this->l1egJetCluster.SetEta(eta);
256  this->l1egJetCluster.SetPhi(phi);
257  this->l1egJetCluster.SetM(mass);
258  }
259  };
260 
261  class simpleL1obj {
262  public:
263  bool stale = false; // Hits become stale once used in clustering algorithm to prevent overlap in clusters
265  false; // L1EGs become associated with a tower to find highest ET total for seeding jets
266  bool passesStandaloneSS = false; // Store whether any of the portions of a WP are passed
267  bool passesStandaloneIso = false; // Store whether any of the portions of a WP are passed
268  bool passesTrkMatchSS = false; // Store whether any of the portions of a WP are passed
269  bool passesTrkMatchIso = false; // Store whether any of the portions of a WP are passed
271 
272  void SetP4(double pt, double eta, double phi, double mass) {
273  this->p4.SetPt(pt);
274  this->p4.SetEta(eta);
275  this->p4.SetPhi(phi);
276  this->p4.SetM(mass);
277  }
278  inline float pt() const { return p4.pt(); };
279  inline float eta() const { return p4.eta(); };
280  inline float phi() const { return p4.phi(); };
281  inline float M() const { return p4.M(); };
282  inline reco::Candidate::PolarLorentzVector GetP4() const { return p4; };
283  };
284 
286  public:
287  int towerIEta = -99;
288  int towerIPhi = -99;
289  float towerEta = -99;
290  float towerPhi = -99;
291  float ecalTowerEt = 0.;
292  float hcalTowerEt = 0.;
293  float l1egTowerEt = 0.;
294  float total_tower_et = 0.;
295  //float total_tower_plus_L1EGs_et=0.;
296  bool stale = false; // Hits become stale once used in clustering algorithm to prevent overlap in clusters
297  bool isBarrel = true; // Defaults to a barrel hit
298 
299  // L1EG info
300  int nL1eg = 0;
301  int l1egTrkSS = 0;
302  int l1egTrkIso = 0;
305  };
306 };
307 
309  : HcalTpEtMin(iConfig.getParameter<double>("HcalTpEtMin")), // Should default to 0 MeV
310  EcalTpEtMin(iConfig.getParameter<double>("EcalTpEtMin")), // Should default to 0 MeV
311  HGCalHadTpEtMin(iConfig.getParameter<double>("HGCalHadTpEtMin")), // Should default to 0 MeV
312  HGCalEmTpEtMin(iConfig.getParameter<double>("HGCalEmTpEtMin")), // Should default to 0 MeV
313  HFTpEtMin(iConfig.getParameter<double>("HFTpEtMin")), // Should default to 0 MeV
314  EtMinForSeedHit(iConfig.getParameter<double>("EtMinForSeedHit")), // Should default to 2.5 GeV
315  EtMinForCollection(iConfig.getParameter<double>("EtMinForCollection")), // Testing 10 GeV
316  EtMinForTauCollection(iConfig.getParameter<double>("EtMinForTauCollection")), // Testing 10 GeV
317  jetPtBins(iConfig.getParameter<std::vector<double>>("jetPtBins")),
318  emFractionBinsBarrel(iConfig.getParameter<std::vector<double>>("emFractionBinsBarrel")),
319  absEtaBinsBarrel(iConfig.getParameter<std::vector<double>>("absEtaBinsBarrel")),
320  jetCalibrationsBarrel(iConfig.getParameter<std::vector<double>>("jetCalibrationsBarrel")),
321  emFractionBinsHGCal(iConfig.getParameter<std::vector<double>>("emFractionBinsHGCal")),
322  absEtaBinsHGCal(iConfig.getParameter<std::vector<double>>("absEtaBinsHGCal")),
323  jetCalibrationsHGCal(iConfig.getParameter<std::vector<double>>("jetCalibrationsHGCal")),
324  emFractionBinsHF(iConfig.getParameter<std::vector<double>>("emFractionBinsHF")),
325  absEtaBinsHF(iConfig.getParameter<std::vector<double>>("absEtaBinsHF")),
326  jetCalibrationsHF(iConfig.getParameter<std::vector<double>>("jetCalibrationsHF")),
327  tauPtBins(iConfig.getParameter<std::vector<double>>("tauPtBins")),
328  tauAbsEtaBinsBarrel(iConfig.getParameter<std::vector<double>>("tauAbsEtaBinsBarrel")),
329  tauCalibrationsBarrel(iConfig.getParameter<std::vector<double>>("tauCalibrationsBarrel")),
330  tauL1egInfoBarrel(iConfig.getParameter<std::vector<edm::ParameterSet>>("tauL1egInfoBarrel")),
331  tauAbsEtaBinsHGCal(iConfig.getParameter<std::vector<double>>("tauAbsEtaBinsHGCal")),
332  tauCalibrationsHGCal(iConfig.getParameter<std::vector<double>>("tauCalibrationsHGCal")),
333  tauL1egInfoHGCal(iConfig.getParameter<std::vector<edm::ParameterSet>>("tauL1egInfoHGCal")),
334  debug(iConfig.getParameter<bool>("debug")),
335  l1TowerToken_(consumes<l1tp2::CaloTowerCollection>(iConfig.getParameter<edm::InputTag>("l1CaloTowers")))
336 
337 {
338  produces<l1tp2::CaloJetsCollection>("L1CaloJetsNoCuts");
339  //produces<l1tp2::CaloJetsCollection>("L1CaloJetsWithCuts");
340  //produces<l1extra::L1JetParticleCollection>("L1CaloClusterCollectionWithCuts");
341  produces<BXVector<l1t::Jet>>("L1CaloJetCollectionBXV");
342  produces<BXVector<l1t::Tau>>("L1CaloTauCollectionBXV");
343 
344  if (debug) {
345  LogDebug("L1CaloJetProducer") << " HcalTpEtMin = " << HcalTpEtMin << " EcalTpEtMin = " << EcalTpEtMin << "\n";
346  }
347 
348  // Fill the calibration 3D vector
349  // Dimension 1 is AbsEta bin
350  // Dimension 2 is EM Fraction bin
351  // Dimension 3 is jet pT bin which is filled with the actual callibration value
352  // size()-1 b/c the inputs have lower and upper bounds
353  // Do Barrel, then HGCal, then HF
354  int index = 0;
355  //calibrations[em_frac][abs_eta].push_back( jetCalibrationsBarrel.at(index) );
356  for (unsigned int abs_eta = 0; abs_eta < absEtaBinsBarrel.size() - 1; abs_eta++) {
357  std::vector<std::vector<double>> em_bins;
358  for (unsigned int em_frac = 0; em_frac < emFractionBinsBarrel.size() - 1; em_frac++) {
359  std::vector<double> pt_bin_calibs;
360  for (unsigned int pt = 0; pt < jetPtBins.size() - 1; pt++) {
361  pt_bin_calibs.push_back(jetCalibrationsBarrel.at(index));
362  index++;
363  }
364  em_bins.push_back(pt_bin_calibs);
365  }
366  calibrationsBarrel.push_back(em_bins);
367  }
368  if (debug) {
369  LogDebug("L1CaloJetProducer") << " Loading Barrel calibrations: Loaded " << index
370  << " values vs. size() of input calibration file: "
371  << int(jetCalibrationsBarrel.size()) << "\n";
372  }
373 
374  index = 0;
375  //calibrations[em_frac][abs_eta].push_back( jetCalibrationsHGCal.at(index) );
376  for (unsigned int abs_eta = 0; abs_eta < absEtaBinsHGCal.size() - 1; abs_eta++) {
377  std::vector<std::vector<double>> em_bins;
378  for (unsigned int em_frac = 0; em_frac < emFractionBinsHGCal.size() - 1; em_frac++) {
379  std::vector<double> pt_bin_calibs;
380  for (unsigned int pt = 0; pt < jetPtBins.size() - 1; pt++) {
381  pt_bin_calibs.push_back(jetCalibrationsHGCal.at(index));
382  index++;
383  }
384  em_bins.push_back(pt_bin_calibs);
385  }
386  calibrationsHGCal.push_back(em_bins);
387  }
388  if (debug) {
389  LogDebug("L1CaloJetProducer") << " Loading HGCal calibrations: Loaded " << index
390  << " values vs. size() of input calibration file: "
391  << int(jetCalibrationsHGCal.size()) << "\n";
392  }
393 
394  index = 0;
395  //calibrations[em_frac][abs_eta].push_back( jetCalibrationsHF.at(index) );
396  for (unsigned int abs_eta = 0; abs_eta < absEtaBinsHF.size() - 1; abs_eta++) {
397  std::vector<std::vector<double>> em_bins;
398  for (unsigned int em_frac = 0; em_frac < emFractionBinsHF.size() - 1; em_frac++) {
399  std::vector<double> pt_bin_calibs;
400  for (unsigned int pt = 0; pt < jetPtBins.size() - 1; pt++) {
401  pt_bin_calibs.push_back(jetCalibrationsHF.at(index));
402  index++;
403  }
404  em_bins.push_back(pt_bin_calibs);
405  }
406  calibrationsHF.push_back(em_bins);
407  }
408  if (debug) {
409  LogDebug("L1CaloJetProducer") << " Loading HF calibrations: Loaded " << index
410  << " values vs. size() of input calibration file: " << int(jetCalibrationsHF.size())
411  << "\n";
412  }
413 
414  // Load Tau L1EG-base calibration info into maps
415  for (auto &first : tauL1egInfoBarrel) {
416  if (debug) {
417  LogDebug("L1CaloJetProducer") << " barrel l1egCount = " << first.getParameter<double>("l1egCount") << "\n";
418  for (auto &em_frac : first.getParameter<std::vector<double>>("l1egEmFractions")) {
419  LogDebug("L1CaloJetProducer") << " - EM = " << em_frac << "\n";
420  }
421  }
422  double l1egCount = first.getParameter<double>("l1egCount");
423  std::vector<double> l1egEmFractions = first.getParameter<std::vector<double>>("l1egEmFractions");
425  tauL1egValuesBarrel.push_back(l1egCount);
426  }
428  for (auto &first : tauL1egInfoHGCal) {
429  if (debug) {
430  LogDebug("L1CaloJetProducer") << " hgcal l1egCount = " << first.getParameter<double>("l1egCount") << "\n";
431  for (auto &em_frac : first.getParameter<std::vector<double>>("l1egEmFractions")) {
432  LogDebug("L1CaloJetProducer") << " - EM = " << em_frac << "\n";
433  }
434  }
435  double l1egCount = first.getParameter<double>("l1egCount");
436  std::vector<double> l1egEmFractions = first.getParameter<std::vector<double>>("l1egEmFractions");
438  tauL1egValuesHGCal.push_back(l1egCount);
439  }
441  // Fill the calibration 4D vector
442  // Dimension 1 is AbsEta bin
443  // Dimension 2 is L1EG count
444  // Dimension 3 is EM Fraction bin
445  // Dimension 4 is tau pT bin which is filled with the actual callibration value
446  // size()-1 b/c the inputs have lower and upper bounds (except L1EG b/c that is a cound)
447  // Do Barrel, then HGCal
448  //
449  // Note to future developers: be very concious of the order in which the calibrations are printed
450  // out in tool which makse the cfg files. You need to match that exactly when loading them and
451  // using the calibrations below.
452  index = 0;
453  for (unsigned int abs_eta = 0; abs_eta < tauAbsEtaBinsBarrel.size() - 1; abs_eta++) {
454  std::vector<std::vector<std::vector<double>>> l1eg_bins;
455  for (auto &l1eg_info : tauL1egInfoMapBarrel) {
456  std::vector<std::vector<double>> em_bins;
457  for (unsigned int em_frac = 0; em_frac < l1eg_info.second.size() - 1; em_frac++) {
458  std::vector<double> pt_bin_calibs;
459  for (unsigned int pt = 0; pt < tauPtBins.size() - 1; pt++) {
460  pt_bin_calibs.push_back(tauCalibrationsBarrel.at(index));
461  index++;
462  }
463  em_bins.push_back(pt_bin_calibs);
464  }
465  l1eg_bins.push_back(em_bins);
466  }
467  tauPtCalibrationsBarrel.push_back(l1eg_bins);
468  }
469  if (debug) {
470  LogDebug("L1CaloJetProducer") << " Loading Barrel calibrations: Loaded " << index
471  << " values vs. size() of input calibration file: "
472  << int(tauCalibrationsBarrel.size()) << "\n";
473  }
474 
475  index = 0;
476  for (unsigned int abs_eta = 0; abs_eta < tauAbsEtaBinsHGCal.size() - 1; abs_eta++) {
477  std::vector<std::vector<std::vector<double>>> l1eg_bins;
478  for (const auto &l1eg_info : tauL1egInfoMapHGCal) {
479  std::vector<std::vector<double>> em_bins;
480  for (unsigned int em_frac = 0; em_frac < l1eg_info.second.size() - 1; em_frac++) {
481  std::vector<double> pt_bin_calibs;
482  for (unsigned int pt = 0; pt < tauPtBins.size() - 1; pt++) {
483  pt_bin_calibs.push_back(tauCalibrationsHGCal.at(index));
484  index++;
485  }
486  em_bins.push_back(pt_bin_calibs);
487  }
488  l1eg_bins.push_back(em_bins);
489  }
490  tauPtCalibrationsHGCal.push_back(l1eg_bins);
491  }
492  if (debug) {
493  LogDebug("L1CaloJetProducer") << " Loading HGCal calibrations: Loaded " << index
494  << " values vs. size() of input calibration file: "
495  << int(tauCalibrationsHGCal.size()) << "\n";
496  }
497 
498  isoTauBarrel.SetParameter(0, 0.30);
499  isoTauBarrel.SetParameter(1, 0.31);
500  isoTauBarrel.SetParameter(2, 0.040);
501  isoTauHGCal.SetParameter(0, 0.34);
502  isoTauHGCal.SetParameter(1, 0.35);
503  isoTauHGCal.SetParameter(2, 0.051);
504 }
505 
507  // Output collections
508  std::unique_ptr<l1tp2::CaloJetsCollection> L1CaloJetsNoCuts(new l1tp2::CaloJetsCollection);
509  //std::unique_ptr<l1tp2::CaloJetsCollection> L1CaloJetsWithCuts( new l1tp2::CaloJetsCollection );
510  //std::unique_ptr<l1extra::L1JetParticleCollection> L1CaloClusterCollectionWithCuts( new l1extra::L1JetParticleCollection );
511  std::unique_ptr<BXVector<l1t::Jet>> L1CaloJetCollectionBXV(new l1t::JetBxCollection);
512  std::unique_ptr<BXVector<l1t::Tau>> L1CaloTauCollectionBXV(new l1t::TauBxCollection);
513 
514  // Load the ECAL+HCAL tower sums coming from L1EGammaCrystalsEmulatorProducer.cc
515  std::vector<SimpleCaloHit> l1CaloTowers;
516 
518  for (auto &hit : *l1CaloTowerHandle.product()) {
519  SimpleCaloHit l1Hit;
520  l1Hit.ecalTowerEt = hit.ecalTowerEt();
521  l1Hit.hcalTowerEt = hit.hcalTowerEt();
522  l1Hit.l1egTowerEt = hit.l1egTowerEt();
523  // Add min ET thresholds for tower ET
524  if (l1Hit.ecalTowerEt < EcalTpEtMin)
525  l1Hit.ecalTowerEt = 0.0;
526  if (l1Hit.hcalTowerEt < HcalTpEtMin)
527  l1Hit.hcalTowerEt = 0.0;
528  l1Hit.total_tower_et = l1Hit.ecalTowerEt + l1Hit.hcalTowerEt + l1Hit.l1egTowerEt;
529  l1Hit.towerIEta = hit.towerIEta();
530  l1Hit.towerIPhi = hit.towerIPhi();
531  l1Hit.nL1eg = hit.nL1eg();
532  l1Hit.l1egTrkSS = hit.l1egTrkSS();
533  l1Hit.l1egTrkIso = hit.l1egTrkIso();
534  l1Hit.l1egStandaloneSS = hit.l1egStandaloneSS();
535  l1Hit.l1egStandaloneIso = hit.l1egStandaloneIso();
536  l1Hit.isBarrel = hit.isBarrel();
537 
538  // FIXME There is an error in the L1EGammaCrystalsEmulatorProducer.cc which is
539  // returning towers with minimal ECAL energy, and no HCAL energy with these
540  // iEta/iPhi coordinates and eta = -88.653152 and phi = -99.000000.
541  // Skip these for the time being until the upstream code has been debugged
542  if ((int)l1Hit.towerIEta == -1016 && (int)l1Hit.towerIPhi == -962)
543  continue;
544 
545  l1Hit.towerEta = hit.towerEta();
546  l1Hit.towerPhi = hit.towerPhi();
547  l1CaloTowers.push_back(l1Hit);
548  if (debug) {
549  LogDebug("L1CaloJetProducer") << " Tower iEta " << (int)l1Hit.towerIEta << " iPhi " << (int)l1Hit.towerIPhi
550  << " eta " << l1Hit.towerEta << " phi " << l1Hit.towerPhi << " ecal_et "
551  << l1Hit.ecalTowerEt << " hacl_et " << l1Hit.hcalTowerEt << " total_et "
552  << l1Hit.total_tower_et << "\n";
553  }
554  }
555 
556  // Sort the ECAL+HCAL+L1EGs tower sums based on total ET
558  return a.total_tower_et > b.total_tower_et;
559  });
560 
561  /**************************************************************************
562  * Begin with making CaloJets in 9x9 grid based on all energy not included in L1EG Objs.
563  * For reference, Run-I used 12x12 grid and Stage-2/Phase-I used 9x9 grid.
564  * We plan to further study this choice and possibly move towards a more circular shape
565  * Create jetCluster within 9x9 of highest ET seed tower.
566  * 9 trigger towers contains all of an ak-0.4 jets, but overshoots on the corners.
567  ******************************************************************************/
568 
569  // Experimental parameters, don't want to bother with hardcoding them in data format
570  std::map<std::string, float> params;
571 
572  std::vector<l1CaloJetObj> l1CaloJetObjs;
573 
574  // Count the number of unused HCAL TPs so we can stop while loop after done.
575  // Clustering can also stop once there are no seed hits >= EtMinForSeedHit
576  int n_towers = l1CaloTowers.size();
577  int n_stale = 0;
578  bool caloJetClusteringFinished = false;
579  while (!caloJetClusteringFinished && n_towers != n_stale) {
580  l1CaloJetObj caloJetObj;
581  caloJetObj.Init();
582 
583  // First find highest ET ECAL+HCAL+L1EGs tower and use to seed the 9x9 Jet
584  int cnt = 0;
585  for (auto &l1CaloTower : l1CaloTowers) {
586  cnt++;
587  if (l1CaloTower.stale)
588  continue; // skip l1CaloTowers which are already used
589 
590  if (caloJetObj.jetClusterET == 0.0) // this is the first l1CaloTower to seed the jet
591  {
592  // Check if the leading unused tower has ET < min for seeding a jet.
593  // If so, stop jet clustering
594  if (l1CaloTower.total_tower_et < EtMinForSeedHit) {
595  caloJetClusteringFinished = true;
596  continue;
597  }
598  l1CaloTower.stale = true;
599  n_stale++;
600 
601  // Set seed location needed for delta iEta/iPhi, eta/phi comparisons later
602  if (l1CaloTower.isBarrel)
603  caloJetObj.barrelSeeded = true;
604  else
605  caloJetObj.barrelSeeded = false;
606 
607  // 3 4-vectors for ECAL, HCAL, ECAL+HCAL for adding together
609  l1CaloTower.hcalTowerEt, l1CaloTower.towerEta, l1CaloTower.towerPhi, 0.);
611  l1CaloTower.ecalTowerEt, l1CaloTower.towerEta, l1CaloTower.towerPhi, 0.);
613  l1CaloTower.l1egTowerEt, l1CaloTower.towerEta, l1CaloTower.towerPhi, 0.);
615  l1CaloTower.total_tower_et, l1CaloTower.towerEta, l1CaloTower.towerPhi, 0.);
616 
617  if (hcalP4.pt() > 0) {
618  caloJetObj.hcal_nHits++;
619  caloJetObj.hcalJetCluster += hcalP4;
620  caloJetObj.hcalJetClusterET += l1CaloTower.hcalTowerEt;
621  }
622  if (ecalP4.pt() > 0) {
623  caloJetObj.ecal_nHits++;
624  caloJetObj.ecalJetCluster += ecalP4;
625  caloJetObj.ecalJetClusterET += l1CaloTower.ecalTowerEt;
626  }
627  if (l1egP4.pt() > 0) {
628  caloJetObj.l1eg_nHits++;
629  caloJetObj.l1egJetCluster += l1egP4;
630  caloJetObj.l1egJetClusterET += l1CaloTower.l1egTowerEt;
631  caloJetObj.l1eg_nL1EGs += l1CaloTower.nL1eg;
632 
633  caloJetObj.l1eg_nL1EGs_standaloneSS += l1CaloTower.l1egStandaloneSS;
634  caloJetObj.l1eg_nL1EGs_standaloneIso += l1CaloTower.l1egStandaloneIso;
635  caloJetObj.l1eg_nL1EGs_trkMatchSS += l1CaloTower.l1egTrkSS;
636  caloJetObj.l1eg_nL1EGs_trkMatchIso += l1CaloTower.l1egTrkIso;
637 
638  if (l1CaloTower.isBarrel) {
639  // For decay mode related checks with CaloTaus
640  // only applicable in the barrel at the moment:
641  // l1eg pt, HCAL ET, ECAL ET, dEta, dPhi, trkSS, trkIso, standaloneSS, standaloneIso
642  std::vector<float> l1EG_info = {float(l1egP4.pt()),
643  float(hcalP4.pt()),
644  float(ecalP4.pt()),
645  0.,
646  0.,
647  float(l1CaloTower.l1egTrkSS),
648  float(l1CaloTower.l1egTrkIso),
649  float(l1CaloTower.l1egStandaloneSS),
650  float(l1CaloTower.l1egStandaloneIso)};
651  if (l1EG_info[1] / (l1EG_info[0] + l1EG_info[2]) < 0.25) {
652  caloJetObj.n_l1eg_HoverE_LessThreshold++;
653  }
654  caloJetObj.associated_l1EGs_.push_back(l1EG_info);
655  }
656  }
657  if (totalP4.pt() > 0) {
658  caloJetObj.total_nHits++;
659  caloJetObj.jetCluster += totalP4;
660  caloJetObj.jetClusterET += l1CaloTower.total_tower_et;
661  caloJetObj.seedTower += totalP4;
662  caloJetObj.seedTowerET += l1CaloTower.total_tower_et;
663  }
664 
665  caloJetObj.seed_iEta = l1CaloTower.towerIEta;
666  caloJetObj.seed_iPhi = l1CaloTower.towerIPhi;
667 
668  if (debug) {
669  LogDebug("L1CaloJetProducer") << " -- hit " << cnt << " , seeding input p4 pt "
670  << l1CaloTower.total_tower_et << " eta " << l1CaloTower.towerEta << " phi "
671  << l1CaloTower.towerPhi << "\n";
672  LogDebug("L1CaloJetProducer") << " -- hit " << cnt << " , seeding input2 p4 pt " << totalP4.pt() << " eta "
673  << totalP4.eta() << " phi " << totalP4.phi() << "\n";
674  LogDebug("L1CaloJetProducer") << " -- hit " << cnt << " seeding reslt tot p4 pt " << caloJetObj.jetClusterET
675  << " eta " << caloJetObj.jetCluster.eta() << " phi "
676  << caloJetObj.jetCluster.phi() << "\n";
677  }
678 
679  // Need to add the seed energy to the dR rings
680  caloJetObj.hcal_seed += hcalP4.pt();
681  //caloJetObj.hcal_3x3 += hcalP4.pt();
682  caloJetObj.hcal_3x5 += hcalP4.pt();
683  //caloJetObj.hcal_5x5 += hcalP4.pt();
684  //caloJetObj.hcal_5x7 += hcalP4.pt();
685  caloJetObj.hcal_7x7 += hcalP4.pt();
686  caloJetObj.ecal_seed += ecalP4.pt();
687  //caloJetObj.ecal_3x3 += ecalP4.pt();
688  caloJetObj.ecal_3x5 += ecalP4.pt();
689  //caloJetObj.ecal_5x5 += ecalP4.pt();
690  //caloJetObj.ecal_5x7 += ecalP4.pt();
691  caloJetObj.ecal_7x7 += ecalP4.pt();
692  caloJetObj.l1eg_seed += l1egP4.pt();
693  //caloJetObj.l1eg_3x3 += l1egP4.pt();
694  caloJetObj.l1eg_3x5 += l1egP4.pt();
695  //caloJetObj.l1eg_5x5 += l1egP4.pt();
696  //caloJetObj.l1eg_5x7 += l1egP4.pt();
697  caloJetObj.l1eg_7x7 += l1egP4.pt();
698  caloJetObj.total_seed += totalP4.pt();
699  //caloJetObj.total_3x3 += totalP4.pt();
700  caloJetObj.total_3x5 += totalP4.pt();
701  //caloJetObj.total_5x5 += totalP4.pt();
702  //caloJetObj.total_5x7 += totalP4.pt();
703  caloJetObj.total_7x7 += totalP4.pt();
704 
705  // Some discrimination vars, 2x2s and 2x3 including central seed
706  //caloJetObj.hcal_2x3 += hcalP4.pt();
707  //caloJetObj.hcal_2x3_1 += hcalP4.pt();
708  //caloJetObj.hcal_2x3_2 += hcalP4.pt();
709  //caloJetObj.ecal_2x3 += ecalP4.pt();
710  //caloJetObj.ecal_2x3_1 += ecalP4.pt();
711  //caloJetObj.ecal_2x3_2 += ecalP4.pt();
712  //caloJetObj.l1eg_2x3 += l1egP4.pt();
713  //caloJetObj.l1eg_2x3_1 += l1egP4.pt();
714  //caloJetObj.l1eg_2x3_2 += l1egP4.pt();
715  //caloJetObj.total_2x3 += totalP4.pt();
716  //caloJetObj.total_2x3_1 += totalP4.pt();
717  //caloJetObj.total_2x3_2 += totalP4.pt();
718 
719  //caloJetObj.hcal_2x2 += hcalP4.pt();
720  //caloJetObj.hcal_2x2_1 += hcalP4.pt();
721  //caloJetObj.hcal_2x2_2 += hcalP4.pt();
722  //caloJetObj.hcal_2x2_3 += hcalP4.pt();
723  //caloJetObj.hcal_2x2_4 += hcalP4.pt();
724  //caloJetObj.ecal_2x2 += ecalP4.pt();
725  //caloJetObj.ecal_2x2_1 += ecalP4.pt();
726  //caloJetObj.ecal_2x2_2 += ecalP4.pt();
727  //caloJetObj.ecal_2x2_3 += ecalP4.pt();
728  //caloJetObj.ecal_2x2_4 += ecalP4.pt();
729  //caloJetObj.l1eg_2x2 += l1egP4.pt();
730  //caloJetObj.l1eg_2x2_1 += l1egP4.pt();
731  //caloJetObj.l1eg_2x2_2 += l1egP4.pt();
732  //caloJetObj.l1eg_2x2_3 += l1egP4.pt();
733  //caloJetObj.l1eg_2x2_4 += l1egP4.pt();
734  //caloJetObj.total_2x2 += totalP4.pt();
735  //caloJetObj.total_2x2_1 += totalP4.pt();
736  //caloJetObj.total_2x2_2 += totalP4.pt();
737  //caloJetObj.total_2x2_3 += totalP4.pt();
738  //caloJetObj.total_2x2_4 += totalP4.pt();
739  continue;
740  }
741 
742  // Unused l1CaloTowers which are not the initial seed
743  // Depending on seed and tower locations calculate iEta/iPhi or eta/phi comparisons.
744  // The defaults of 99 will automatically fail comparisons for the incorrect regions.
745  int hit_iPhi = 99;
746  int d_iEta = 99;
747  int d_iPhi = 99;
748  float d_eta = 99;
749  float d_phi = 99;
750  if (caloJetObj.barrelSeeded && l1CaloTower.isBarrel) // use iEta/iPhi comparisons
751  {
752  hit_iPhi = l1CaloTower.towerIPhi;
753  d_iEta = tower_diEta(caloJetObj.seed_iEta, l1CaloTower.towerIEta);
754  d_iPhi = tower_diPhi(caloJetObj.seed_iPhi, hit_iPhi);
755  } else // either seed or tower are in HGCal or HF, use eta/phi
756  {
757  d_eta = caloJetObj.seedTower.eta() - l1CaloTower.towerEta;
758  d_phi = reco::deltaPhi(caloJetObj.seedTower.phi(), l1CaloTower.towerPhi);
759  }
760 
761  // 7x7 HCAL Trigger Towers
762  // If seeded in barrel and hit is barrel then we can compare iEta/iPhi, else need to use eta/phi
763  // in HGCal / transition region
764  if ((abs(d_iEta) <= 3 && abs(d_iPhi) <= 3) || (std::abs(d_eta) < 0.3 && std::abs(d_phi) < 0.3)) {
765  l1CaloTower.stale = true;
766  n_stale++;
767 
768  // 3 4-vectors for ECAL, HCAL, ECAL+HCAL for adding together
770  l1CaloTower.hcalTowerEt, l1CaloTower.towerEta, l1CaloTower.towerPhi, 0.);
772  l1CaloTower.ecalTowerEt, l1CaloTower.towerEta, l1CaloTower.towerPhi, 0.);
774  l1CaloTower.l1egTowerEt, l1CaloTower.towerEta, l1CaloTower.towerPhi, 0.);
776  l1CaloTower.total_tower_et, l1CaloTower.towerEta, l1CaloTower.towerPhi, 0.);
777 
778  if (hcalP4.pt() > 0) {
779  caloJetObj.hcal_nHits++;
780  caloJetObj.hcalJetCluster += hcalP4;
781  caloJetObj.hcalJetClusterET += l1CaloTower.hcalTowerEt;
782  }
783  if (ecalP4.pt() > 0) {
784  caloJetObj.ecal_nHits++;
785  caloJetObj.ecalJetCluster += ecalP4;
786  caloJetObj.ecalJetClusterET += l1CaloTower.ecalTowerEt;
787  }
788  if (l1egP4.pt() > 0) {
789  caloJetObj.l1eg_nHits++;
790  caloJetObj.l1egJetCluster += l1egP4;
791  caloJetObj.l1egJetClusterET += l1CaloTower.l1egTowerEt;
792  caloJetObj.l1eg_nL1EGs += l1CaloTower.nL1eg;
793  }
794  if (totalP4.pt() > 0) {
795  caloJetObj.total_nHits++;
796  caloJetObj.jetCluster += totalP4;
797  caloJetObj.jetClusterET += l1CaloTower.total_tower_et;
798  }
799 
800  if (debug) {
801  LogDebug("L1CaloJetProducer") << " ---- hit " << cnt << " input p4 pt " << totalP4.pt() << " eta "
802  << totalP4.eta() << " phi " << totalP4.phi() << "\n";
803  LogDebug("L1CaloJetProducer") << " ---- hit " << cnt << " resulting p4 pt " << caloJetObj.jetClusterET
804  << " eta " << caloJetObj.jetCluster.eta() << " phi "
805  << caloJetObj.jetCluster.phi() << "\n";
806  }
807 
808  if ((abs(d_iEta) == 0 && abs(d_iPhi) == 0) || (std::abs(d_eta) < 0.043 && std::abs(d_phi) < 0.043)) {
809  caloJetObj.hcal_seed += hcalP4.pt();
810  caloJetObj.ecal_seed += ecalP4.pt();
811  caloJetObj.l1eg_seed += l1egP4.pt();
812  caloJetObj.total_seed += totalP4.pt();
813  }
814  //if ( (abs( d_iEta ) <= 1 && abs( d_iPhi ) <= 1) ||
815  // ( std::abs( d_eta ) < 0.13 && std::abs( d_phi ) < 0.13 ) )
816  //{
817  // caloJetObj.hcal_3x3 += hcalP4.pt();
818  // caloJetObj.ecal_3x3 += ecalP4.pt();
819  // caloJetObj.l1eg_3x3 += l1egP4.pt();
820  // caloJetObj.total_3x3 += totalP4.pt();
821  //}
822  if ((abs(d_iEta) <= 1 && abs(d_iPhi) <= 2) || (std::abs(d_eta) < 0.13 && std::abs(d_phi) < 0.22)) {
823  caloJetObj.hcal_3x5 += hcalP4.pt();
824  caloJetObj.ecal_3x5 += ecalP4.pt();
825  caloJetObj.l1eg_3x5 += l1egP4.pt();
826  caloJetObj.total_3x5 += totalP4.pt();
827 
828  // Do this for 3x5 only
829  if (l1egP4.pt() > 0) {
830  caloJetObj.l1eg_nL1EGs_standaloneSS += l1CaloTower.l1egStandaloneSS;
831  caloJetObj.l1eg_nL1EGs_standaloneIso += l1CaloTower.l1egStandaloneIso;
832  caloJetObj.l1eg_nL1EGs_trkMatchSS += l1CaloTower.l1egTrkSS;
833  caloJetObj.l1eg_nL1EGs_trkMatchIso += l1CaloTower.l1egTrkIso;
834 
835  // For decay mode related checks with CaloTaus
836  // only applicable in the barrel at the moment:
837  // l1eg pt, HCAL ET, ECAL ET, d_iEta, d_iPhi, trkSS, trkIso, standaloneSS, standaloneIso
838  std::vector<float> l1EG_info = {float(l1egP4.pt()),
839  float(hcalP4.pt()),
840  float(ecalP4.pt()),
841  float(d_iEta),
842  float(d_iPhi),
843  float(l1CaloTower.l1egTrkSS),
844  float(l1CaloTower.l1egTrkIso),
845  float(l1CaloTower.l1egStandaloneSS),
846  float(l1CaloTower.l1egStandaloneIso)};
847  if (l1EG_info[1] / (l1EG_info[0] + l1EG_info[2]) < 0.25) {
848  caloJetObj.n_l1eg_HoverE_LessThreshold++;
849  }
850  caloJetObj.associated_l1EGs_.push_back(l1EG_info);
851  }
852  }
853  //if ( ( abs( d_iEta ) <= 2 && abs( d_iPhi ) <= 2) ||
854  // ( std::abs( d_eta ) < 0.22 && std::abs( d_phi ) < 0.22 ) )
855  //{
856  // caloJetObj.hcal_5x5 += hcalP4.pt();
857  // caloJetObj.ecal_5x5 += ecalP4.pt();
858  // caloJetObj.l1eg_5x5 += l1egP4.pt();
859  // caloJetObj.total_5x5 += totalP4.pt();
860  //}
861  //if ( ( abs( d_iEta ) <= 2 && abs( d_iPhi ) <= 3) ||
862  // ( std::abs( d_eta ) < 0.22 && std::abs( d_phi ) < 0.3 ) )
863  //{
864  // caloJetObj.hcal_5x7 += hcalP4.pt();
865  // caloJetObj.ecal_5x7 += ecalP4.pt();
866  // caloJetObj.l1eg_5x7 += l1egP4.pt();
867  // caloJetObj.total_5x7 += totalP4.pt();
868  //}
869  if ((abs(d_iEta) <= 3 && abs(d_iPhi) <= 3) || (std::abs(d_eta) < 0.3 && std::abs(d_phi) < 0.3)) {
870  caloJetObj.hcal_7x7 += hcalP4.pt();
871  caloJetObj.ecal_7x7 += ecalP4.pt();
872  caloJetObj.l1eg_7x7 += l1egP4.pt();
873  caloJetObj.total_7x7 += totalP4.pt();
874  }
875 
878  //if ( ( d_iEta == 0 || d_iEta == 1 ) && abs(d_iPhi) <= 1 )
879  //{
880  // caloJetObj.hcal_2x3_1 += hcalP4.pt();
881  // caloJetObj.ecal_2x3_1 += ecalP4.pt();
882  // caloJetObj.l1eg_2x3_1 += l1egP4.pt();
883  // caloJetObj.total_2x3_1 += totalP4.pt();
884  //}
885  //if ( ( d_iEta == 0 || d_iEta == -1 ) && abs(d_iPhi) <= 1 )
886  //{
887  // caloJetObj.hcal_2x3_2 += hcalP4.pt();
888  // caloJetObj.ecal_2x3_2 += ecalP4.pt();
889  // caloJetObj.l1eg_2x3_2 += l1egP4.pt();
890  // caloJetObj.total_2x3_2 += totalP4.pt();
891  //}
893  //if ( std::abs( d_eta ) < 0.087 && std::abs( d_phi ) < 0.13 )
894  //{
895  // caloJetObj.hcal_2x3 += hcalP4.pt();
896  // caloJetObj.ecal_2x3 += ecalP4.pt();
897  // caloJetObj.l1eg_2x3 += l1egP4.pt();
898  // caloJetObj.total_2x3 += totalP4.pt();
899  //}
900 
902  //if ( ( d_iEta == 0 || d_iEta == 1 ) && ( d_iPhi == 0 || d_iPhi == 1 ) )
903  //{
904  // caloJetObj.hcal_2x2_1 += hcalP4.pt();
905  // caloJetObj.ecal_2x2_1 += ecalP4.pt();
906  // caloJetObj.l1eg_2x2_1 += l1egP4.pt();
907  // caloJetObj.total_2x2_1 += totalP4.pt();
908  //}
909  //if ( ( d_iEta == 0 || d_iEta == 1 ) && ( d_iPhi == 0 || d_iPhi == -1 ) )
910  //{
911  // caloJetObj.hcal_2x2_2 += hcalP4.pt();
912  // caloJetObj.ecal_2x2_2 += ecalP4.pt();
913  // caloJetObj.l1eg_2x2_2 += l1egP4.pt();
914  // caloJetObj.total_2x2_2 += totalP4.pt();
915  //}
916  //if ( ( d_iEta == 0 || d_iEta == -1 ) && ( d_iPhi == 0 || d_iPhi == 1 ) )
917  //{
918  // caloJetObj.hcal_2x2_3 += hcalP4.pt();
919  // caloJetObj.ecal_2x2_3 += ecalP4.pt();
920  // caloJetObj.l1eg_2x2_3 += l1egP4.pt();
921  // caloJetObj.total_2x2_3 += totalP4.pt();
922  //}
923  //if ( ( d_iEta == 0 || d_iEta == -1 ) && ( d_iPhi == 0 || d_iPhi == -1 ) )
924  //{
925  // caloJetObj.hcal_2x2_4 += hcalP4.pt();
926  // caloJetObj.ecal_2x2_4 += ecalP4.pt();
927  // caloJetObj.l1eg_2x2_4 += l1egP4.pt();
928  // caloJetObj.total_2x2_4 += totalP4.pt();
929  //}
931  //if ( std::abs( d_eta ) < 0.087 && std::abs( d_phi ) < 0.087 )
932  //{
933  // caloJetObj.hcal_2x2 += hcalP4.pt();
934  // caloJetObj.ecal_2x2 += ecalP4.pt();
935  // caloJetObj.l1eg_2x2 += l1egP4.pt();
936  // caloJetObj.total_2x2 += totalP4.pt();
937  //}
938  }
939  }
940 
941  if (caloJetObj.jetClusterET > 0.0) {
942  l1CaloJetObjs.push_back(caloJetObj);
943  }
944 
945  } // end while loop of HCAL TP clustering
946 
947  // Sort JetClusters so we can begin with the highest pt for next step of jet clustering
948  std::sort(begin(l1CaloJetObjs), end(l1CaloJetObjs), [](const l1CaloJetObj &a, const l1CaloJetObj &b) {
949  return a.jetClusterET > b.jetClusterET;
950  });
951 
952  /**************************************************************************
953  * Progress to adding L1EGs built from ECAL TPs 9x9 grid.
954  * Recall, for 9x9 trigger towers gives diameter 0.78
955  ******************************************************************************/
956 
957  // Cluster together the L1EGs around existing HCAL Jet
958  // Cluster within dEta/dPhi 0.4 which is very close to 0.39 = 9x9/2
959  //std::cout << " - Input L1EGs: " << crystalClustersVect.size() << std::endl;
960  for (auto &caloJetObj : l1CaloJetObjs) {
961  params["seed_pt"] = caloJetObj.seedTowerET;
962  params["seed_eta"] = caloJetObj.seedTower.eta();
963  params["seed_phi"] = caloJetObj.seedTower.phi();
964  params["seed_iEta"] = caloJetObj.seed_iEta;
965  params["seed_iPhi"] = caloJetObj.seed_iPhi;
966  params["seed_energy"] = caloJetObj.seedTower.energy();
967 
968  params["hcal_pt"] = caloJetObj.hcalJetClusterET;
969  params["hcal_seed"] = caloJetObj.hcal_seed;
970  //params["hcal_3x3"] = caloJetObj.hcal_3x3;
971  params["hcal_3x5"] = caloJetObj.hcal_3x5;
972  //params["hcal_5x5"] = caloJetObj.hcal_5x5;
973  //params["hcal_5x7"] = caloJetObj.hcal_5x7;
974  params["hcal_7x7"] = caloJetObj.hcal_7x7;
975  //params["hcal_2x3"] = std::max( caloJetObj.hcal_2x3, std::max( caloJetObj.hcal_2x3_1, caloJetObj.hcal_2x3_2 ));
976  //params["hcal_2x2"] = std::max( caloJetObj.hcal_2x2, std::max( caloJetObj.hcal_2x2_1, std::max( caloJetObj.hcal_2x2_2, std::max( caloJetObj.hcal_2x2_3, caloJetObj.hcal_2x2_4 ))));
977  params["hcal_nHits"] = caloJetObj.hcal_nHits;
978 
979  params["ecal_pt"] = caloJetObj.ecalJetClusterET;
980  params["ecal_seed"] = caloJetObj.ecal_seed;
981  //params["ecal_3x3"] = caloJetObj.ecal_3x3;
982  params["ecal_3x5"] = caloJetObj.ecal_3x5;
983  //params["ecal_5x5"] = caloJetObj.ecal_5x5;
984  //params["ecal_5x7"] = caloJetObj.ecal_5x7;
985  params["ecal_7x7"] = caloJetObj.ecal_7x7;
986  //params["ecal_2x3"] = std::max( caloJetObj.ecal_2x3, std::max( caloJetObj.ecal_2x3_1, caloJetObj.ecal_2x3_2 ));
987  //params["ecal_2x2"] = std::max( caloJetObj.ecal_2x2, std::max( caloJetObj.ecal_2x2_1, std::max( caloJetObj.ecal_2x2_2, std::max( caloJetObj.ecal_2x2_3, caloJetObj.ecal_2x2_4 ))));
988  params["ecal_nHits"] = caloJetObj.ecal_nHits;
989 
990  params["l1eg_pt"] = caloJetObj.l1egJetClusterET;
991  params["l1eg_seed"] = caloJetObj.l1eg_seed;
992  //params["l1eg_3x3"] = caloJetObj.l1eg_3x3;
993  params["l1eg_3x5"] = caloJetObj.l1eg_3x5;
994  //params["l1eg_5x5"] = caloJetObj.l1eg_5x5;
995  //params["l1eg_5x7"] = caloJetObj.l1eg_5x7;
996  params["l1eg_7x7"] = caloJetObj.l1eg_7x7;
997  //params["l1eg_2x3"] = std::max( caloJetObj.l1eg_2x3, std::max( caloJetObj.l1eg_2x3_1, caloJetObj.l1eg_2x3_2 ));
998  //params["l1eg_2x2"] = std::max( caloJetObj.l1eg_2x2, std::max( caloJetObj.l1eg_2x2_1, std::max( caloJetObj.l1eg_2x2_2, std::max( caloJetObj.l1eg_2x2_3, caloJetObj.l1eg_2x2_4 ))));
999  params["l1eg_nHits"] = caloJetObj.l1eg_nHits;
1000  params["l1eg_nL1EGs"] = caloJetObj.l1eg_nL1EGs;
1001  params["l1eg_nL1EGs_standaloneSS"] = caloJetObj.l1eg_nL1EGs_standaloneSS;
1002  params["l1eg_nL1EGs_standaloneIso"] = caloJetObj.l1eg_nL1EGs_standaloneIso;
1003  params["l1eg_nL1EGs_trkMatchSS"] = caloJetObj.l1eg_nL1EGs_trkMatchSS;
1004  params["l1eg_nL1EGs_trkMatchIso"] = caloJetObj.l1eg_nL1EGs_trkMatchIso;
1005 
1006  params["total_et"] = caloJetObj.jetClusterET;
1007  params["total_seed"] = caloJetObj.total_seed;
1008  //params["total_3x3"] = caloJetObj.total_3x3;
1009  params["total_3x5"] = caloJetObj.total_3x5;
1010  //params["total_5x5"] = caloJetObj.total_5x5;
1011  //params["total_5x7"] = caloJetObj.total_5x7;
1012  params["total_7x7"] = caloJetObj.total_7x7;
1013  //params["total_2x3"] = std::max( caloJetObj.total_2x3, std::max( caloJetObj.total_2x3_1, caloJetObj.total_2x3_2 ));
1014  //params["total_2x2"] = std::max( caloJetObj.total_2x2, std::max( caloJetObj.total_2x2_1, std::max( caloJetObj.total_2x2_2, std::max( caloJetObj.total_2x2_3, caloJetObj.total_2x2_4 ))));
1015  params["total_nHits"] = caloJetObj.total_nHits;
1016  //params["total_nTowers"] = total_nTowers;
1017 
1019  float hovere = -9;
1020  if (caloJetObj.ecalJetClusterET > 0.0) {
1021  hovere = caloJetObj.hcalJetClusterET / (caloJetObj.ecalJetClusterET + caloJetObj.l1egJetClusterET);
1022  }
1023 
1024  params["jet_pt"] = caloJetObj.jetClusterET;
1025  params["jet_eta"] = caloJetObj.jetCluster.eta();
1026  params["jet_phi"] = caloJetObj.jetCluster.phi();
1027  params["jet_mass"] = caloJetObj.jetCluster.mass();
1028  params["jet_energy"] = caloJetObj.jetCluster.energy();
1029 
1030  // Calibrations
1031  params["hcal_calibration"] =
1032  get_hcal_calibration(params["jet_pt"], params["ecal_pt"], params["l1eg_pt"], params["jet_eta"]);
1033  params["hcal_pt_calibration"] = params["hcal_pt"] * params["hcal_calibration"];
1034  params["jet_pt_calibration"] = params["hcal_pt_calibration"] + params["ecal_pt"] + params["l1eg_pt"];
1035 
1036  // Tau Vars
1037  // The tau pT calibration is applied as a SF to the total raw pT
1038  // in contrast to the jet calibration above
1039  params["tau_pt_calibration_value"] = get_tau_pt_calibration(params["total_3x5"],
1040  params["ecal_3x5"],
1041  params["l1eg_3x5"],
1042  caloJetObj.n_l1eg_HoverE_LessThreshold,
1043  params["jet_eta"]);
1044  params["tau_pt"] = params["total_3x5"] * params["tau_pt_calibration_value"];
1045  params["n_l1eg_HoverE_LessThreshold"] = caloJetObj.n_l1eg_HoverE_LessThreshold;
1046  // Currently, applying the tau_pt calibration to the isolation region as well...
1047  // One could switch to using the calibrated jet_pt instead for the iso region...
1048  // This should be revisited - FIXME?
1049  params["tau_total_iso_et"] = params["jet_pt"] * params["tau_pt_calibration_value"];
1050  params["tau_iso_et"] = (params["jet_pt"] * params["tau_pt_calibration_value"]) - params["tau_pt"];
1051  params["loose_iso_tau_wp"] = float(loose_iso_tau_wp(params["tau_pt"], params["tau_iso_et"], params["jet_eta"]));
1052 
1053  float calibratedPt = -1;
1054  float ECalIsolation = -1; // Need to loop over 7x7 crystals of unclustered energy
1055  float totalPtPUcorr = -1;
1056  l1tp2::CaloJet caloJet(caloJetObj.jetCluster, calibratedPt, hovere, ECalIsolation, totalPtPUcorr);
1057  caloJet.setExperimentalParams(params);
1058  caloJet.setAssociated_l1EGs(caloJetObj.associated_l1EGs_);
1059 
1060  // Only store jets passing ET threshold
1061  if (params["jet_pt_calibration"] >= EtMinForCollection) {
1062  L1CaloJetsNoCuts->push_back(caloJet);
1063  //L1CaloJetsWithCuts->push_back( caloJet );
1065  params["jet_pt_calibration"], caloJet.p4().eta(), caloJet.p4().phi(), caloJet.p4().M());
1066  L1CaloJetCollectionBXV->push_back(0, l1t::Jet(jet_p4));
1067 
1068  if (debug)
1069  LogDebug("L1CaloJetProducer") << " Made a Jet, eta " << caloJetObj.jetCluster.eta() << " phi "
1070  << caloJetObj.jetCluster.phi() << " pt " << caloJetObj.jetClusterET
1071  << " calibrated pt " << params["jet_pt_calibration"] << "\n";
1072  }
1073 
1074  // Only store taus passing ET threshold
1075  if (params["tau_pt"] >= EtMinForTauCollection) {
1076  short int tau_ieta = caloJetObj.seed_iEta;
1077  short int tau_iphi = caloJetObj.seed_iPhi;
1078  short int raw_et = params["total_3x5"];
1079  short int iso_et = params["tau_iso_et"];
1080  bool hasEM = false;
1081  if (params["l1eg_3x5"] > 0. || params["ecal_3x5"] > 0.) {
1082  hasEM = true;
1083  }
1084  int tau_qual = int(params["loose_iso_tau_wp"]);
1085 
1087  params["tau_pt"], caloJet.p4().eta(), caloJet.p4().phi(), caloJet.p4().M());
1088  l1t::Tau l1Tau = l1t::Tau(tau_p4, params["tau_pt"], caloJet.p4().eta(), caloJet.p4().phi(), tau_qual, iso_et);
1089  l1Tau.setTowerIEta(tau_ieta);
1090  l1Tau.setTowerIPhi(tau_iphi);
1091  l1Tau.setRawEt(raw_et);
1092  l1Tau.setIsoEt(iso_et);
1093  l1Tau.setHasEM(hasEM);
1094  l1Tau.setIsMerged(false);
1095  L1CaloTauCollectionBXV->push_back(0, l1Tau);
1096 
1097  if (debug) {
1098  LogDebug("L1CaloJetProducer") << " Made a Jet, eta " << l1Tau.eta() << " phi " << l1Tau.phi() << " pt "
1099  << l1Tau.rawEt() << " calibrated pt " << l1Tau.pt() << "\n";
1100  }
1101  }
1102 
1103  } // end jetClusters loop
1104 
1105  iEvent.put(std::move(L1CaloJetsNoCuts), "L1CaloJetsNoCuts");
1106  //iEvent.put(std::move(L1CaloJetsWithCuts), "L1CaloJetsWithCuts" );
1107  //iEvent.put(std::move(L1CaloClusterCollectionWithCuts), "L1CaloClusterCollectionWithCuts" );
1108  iEvent.put(std::move(L1CaloJetCollectionBXV), "L1CaloJetCollectionBXV");
1109  iEvent.put(std::move(L1CaloTauCollectionBXV), "L1CaloTauCollectionBXV");
1110 }
1111 
1112 int L1CaloJetProducer::tower_diPhi(int &iPhi_1, int &iPhi_2) const {
1113  // 360 Crystals in full, 72 towers, half way is 36
1114  int PI = 36;
1115  int result = iPhi_1 - iPhi_2;
1116  while (result > PI)
1117  result -= 2 * PI;
1118  while (result <= -PI)
1119  result += 2 * PI;
1120  return result;
1121 }
1122 
1123 // Added b/c of the iEta jump from +1 to -1 across the barrel mid point
1124 int L1CaloJetProducer::tower_diEta(int &iEta_1, int &iEta_2) const {
1125  // On same side of barrel
1126  if (iEta_1 * iEta_2 > 0)
1127  return iEta_1 - iEta_2;
1128  else
1129  return iEta_1 - iEta_2 - 1;
1130 }
1131 
1134  // Check that pt is > 0 for both or else reco::deltaR returns bogus values
1135  if (p4_1.pt() > 0 && p4_2.pt() > 0) {
1136  return reco::deltaR(p4_1, p4_2);
1137  } else
1138  return -1;
1139 }
1140 
1141 // Apply calibrations to HCAL energy based on EM Fraction, Jet Eta, Jet pT
1143  float &ecal_pt,
1144  float &ecal_L1EG_jet_pt,
1145  float &jet_eta) const {
1146  float em_frac = (ecal_L1EG_jet_pt + ecal_pt) / jet_pt;
1147  float abs_eta = std::abs(jet_eta);
1148  float tmp_jet_pt = jet_pt;
1149  if (tmp_jet_pt > 499)
1150  tmp_jet_pt = 499;
1151 
1152  // Different indices sizes in different calo regions.
1153  // Barrel...
1154  size_t em_index = 0;
1155  size_t eta_index = 0;
1156  size_t pt_index = 0;
1157  float calib = 1.0;
1158  if (abs_eta <= 1.5) {
1159  // Start loop checking 2nd value
1160  for (unsigned int i = 1; i < emFractionBinsBarrel.size(); i++) {
1161  if (em_frac <= emFractionBinsBarrel.at(i))
1162  break;
1163  em_index++;
1164  }
1165 
1166  // Start loop checking 2nd value
1167  for (unsigned int i = 1; i < absEtaBinsBarrel.size(); i++) {
1168  if (abs_eta <= absEtaBinsBarrel.at(i))
1169  break;
1170  eta_index++;
1171  }
1172 
1173  // Start loop checking 2nd value
1174  for (unsigned int i = 1; i < jetPtBins.size(); i++) {
1175  if (tmp_jet_pt <= jetPtBins.at(i))
1176  break;
1177  pt_index++;
1178  }
1179  calib = calibrationsBarrel[eta_index][em_index][pt_index];
1180  } // end Barrel
1181  else if (abs_eta <= 3.0) // HGCal
1182  {
1183  // Start loop checking 2nd value
1184  for (unsigned int i = 1; i < emFractionBinsHGCal.size(); i++) {
1185  if (em_frac <= emFractionBinsHGCal.at(i))
1186  break;
1187  em_index++;
1188  }
1189 
1190  // Start loop checking 2nd value
1191  for (unsigned int i = 1; i < absEtaBinsHGCal.size(); i++) {
1192  if (abs_eta <= absEtaBinsHGCal.at(i))
1193  break;
1194  eta_index++;
1195  }
1196 
1197  // Start loop checking 2nd value
1198  for (unsigned int i = 1; i < jetPtBins.size(); i++) {
1199  if (tmp_jet_pt <= jetPtBins.at(i))
1200  break;
1201  pt_index++;
1202  }
1203  calib = calibrationsHGCal[eta_index][em_index][pt_index];
1204  } // end HGCal
1205  else // HF
1206  {
1207  // Start loop checking 2nd value
1208  for (unsigned int i = 1; i < emFractionBinsHF.size(); i++) {
1209  if (em_frac <= emFractionBinsHF.at(i))
1210  break;
1211  em_index++;
1212  }
1213 
1214  // Start loop checking 2nd value
1215  for (unsigned int i = 1; i < absEtaBinsHF.size(); i++) {
1216  if (abs_eta <= absEtaBinsHF.at(i))
1217  break;
1218  eta_index++;
1219  }
1220 
1221  // Start loop checking 2nd value
1222  for (unsigned int i = 1; i < jetPtBins.size(); i++) {
1223  if (tmp_jet_pt <= jetPtBins.at(i))
1224  break;
1225  pt_index++;
1226  }
1227  calib = calibrationsHF[eta_index][em_index][pt_index];
1228  } // end HF
1229 
1230  return calib;
1231 }
1232 
1233 // Apply calibrations to tau pT based on L1EG info, EM Fraction, Tau Eta, Tau pT
1235  float &tau_pt, float &ecal_pt, float &l1EG_pt, float &n_L1EGs, float &tau_eta) const {
1236  float em_frac = (l1EG_pt + ecal_pt) / tau_pt;
1237  float abs_eta = std::abs(tau_eta);
1238  float tmp_tau_pt = tau_pt;
1239  if (tmp_tau_pt > 199)
1240  tmp_tau_pt = 199;
1241 
1242  // Different indices sizes in different calo regions.
1243  // Barrel...
1244  size_t em_index = 0;
1245  size_t eta_index = 0;
1246  size_t n_L1EG_index = 0;
1247  size_t pt_index = 0;
1248  float calib = 1.0;
1249  // HERE
1250  if (abs_eta <= 1.5) {
1251  // Start loop checking 1st value
1252  for (unsigned int i = 0; i < tauL1egValuesBarrel.size(); i++) {
1253  if (n_L1EGs == tauL1egValuesBarrel.at(i))
1254  break;
1255  if (tauL1egValuesBarrel.at(i) == tauL1egValuesBarrel.back())
1256  break; // to preven incrementing on last one
1257  n_L1EG_index++;
1258  }
1259 
1260  // Find key value pair matching n L1EGs
1261  for (auto &l1eg_info : tauL1egInfoMapBarrel) {
1262  if (l1eg_info.first != double(n_L1EG_index))
1263  continue;
1264  // Start loop checking 2nd value
1265  for (unsigned int i = 1; i < l1eg_info.second.size(); i++) {
1266  if (em_frac <= l1eg_info.second.at(i))
1267  break;
1268  em_index++;
1269  }
1270  }
1271 
1272  // Start loop checking 2nd value
1273  for (unsigned int i = 1; i < tauAbsEtaBinsBarrel.size(); i++) {
1274  if (abs_eta <= tauAbsEtaBinsBarrel.at(i))
1275  break;
1276  eta_index++;
1277  }
1278 
1279  // Start loop checking 2nd value
1280  for (unsigned int i = 1; i < tauPtBins.size(); i++) {
1281  if (tmp_tau_pt <= tauPtBins.at(i))
1282  break;
1283  pt_index++;
1284  }
1285  calib = tauPtCalibrationsBarrel[eta_index][n_L1EG_index][em_index][pt_index];
1286  } // end Barrel
1287  else if (abs_eta <= 3.0) // HGCal
1288  {
1289  // Start loop checking 1st value
1290  for (unsigned int i = 0; i < tauL1egValuesHGCal.size(); i++) {
1291  if (n_L1EGs == tauL1egValuesHGCal.at(i))
1292  break;
1293  if (tauL1egValuesHGCal.at(i) == tauL1egValuesHGCal.back())
1294  break; // to preven incrementing on last one
1295  n_L1EG_index++;
1296  }
1297 
1298  // Find key value pair matching n L1EGs
1299  for (auto &l1eg_info : tauL1egInfoMapHGCal) {
1300  if (l1eg_info.first != double(n_L1EG_index))
1301  continue;
1302  // Start loop checking 2nd value
1303  for (unsigned int i = 1; i < l1eg_info.second.size(); i++) {
1304  if (em_frac <= l1eg_info.second.at(i))
1305  break;
1306  em_index++;
1307  }
1308  }
1309 
1310  // Start loop checking 2nd value
1311  for (unsigned int i = 1; i < tauAbsEtaBinsHGCal.size(); i++) {
1312  if (abs_eta <= tauAbsEtaBinsHGCal.at(i))
1313  break;
1314  eta_index++;
1315  }
1316 
1317  // Start loop checking 2nd value
1318  for (unsigned int i = 1; i < tauPtBins.size(); i++) {
1319  if (tmp_tau_pt <= tauPtBins.at(i))
1320  break;
1321  pt_index++;
1322  }
1323  calib = tauPtCalibrationsHGCal[eta_index][n_L1EG_index][em_index][pt_index];
1324  } // end HGCal
1325  else
1326  return calib;
1327 
1328  return calib;
1329 }
1330 
1331 // Loose IsoTau WP
1332 int L1CaloJetProducer::loose_iso_tau_wp(float &tau_pt, float &tau_iso_et, float &tau_eta) const {
1333  // Fully relaxed above 100 GeV pT
1334  if (tau_pt > 100) {
1335  return 1;
1336  }
1337  // Split by barrel and HGCal
1338  // with Barrel first
1339  if (std::abs(tau_eta) < 1.5) {
1340  if (isoTauBarrel.Eval(tau_pt) >= (tau_iso_et / tau_pt)) {
1341  return 1;
1342  } else {
1343  return 0;
1344  }
1345  }
1346  // HGCal
1347  if (std::abs(tau_eta) < 3.0) {
1348  if (isoTauHGCal.Eval(tau_pt) >= (tau_iso_et / tau_pt)) {
1349  return 1;
1350  } else {
1351  return 0;
1352  }
1353  }
1354  // Beyond HGCal
1355  return 0;
1356 }
1357 
void SetEcalJetClusterP4(double pt, double eta, double phi, double mass)
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
void SetL1EGJetP4(double pt, double eta, double phi, double mass)
int loose_iso_tau_wp(float &tau_pt, float &tau_iso_et, float &tau_eta) const
reco::Candidate::PolarLorentzVector hcalJetCluster
std::vector< double > tauAbsEtaBinsBarrel
std::map< double, std::vector< double > > tauL1egInfoMapBarrel
void setExperimentalParams(const std::map< std::string, float > &params)
Definition: CaloJet.h:25
std::vector< double > absEtaBinsHF
std::vector< double > jetCalibrationsHF
Definition: Tau.h:20
std::vector< double > jetCalibrationsBarrel
void setAssociated_l1EGs(const std::vector< std::vector< float >> l1EGs)
Definition: CaloJet.h:26
T const * product() const
Definition: Handle.h:70
std::vector< std::vector< float > > associated_l1EGs_
int tower_diPhi(int &iPhi_1, int &iPhi_2) const
edm::EDGetTokenT< l1tp2::CaloTowerCollection > l1TowerToken_
void SetSeedP4(double pt, double eta, double phi, double mass)
float get_deltaR(reco::Candidate::PolarLorentzVector &p4_1, reco::Candidate::PolarLorentzVector &p4_2) const
float get_hcal_calibration(float &jet_pt, float &ecal_pt, float &ecal_L1EG_jet_pt, float &jet_eta) const
std::vector< double > emFractionBinsHGCal
float get_tau_pt_calibration(float &tau_pt, float &ecal_pt, float &l1EG_pt, float &n_L1EGs, float &tau_eta) const
std::vector< std::vector< std::vector< double > > > calibrationsBarrel
const LorentzVector & p4() const final
four-momentum Lorentz vector
std::vector< double > tauAbsEtaBinsHGCal
Definition: Jet.h:20
void SetP4(double pt, double eta, double phi, double mass)
int iEvent
Definition: GenABIO.cc:224
std::vector< l1tp2::CaloJet > CaloJetsCollection
Definition: CaloJet.h:57
std::vector< double > emFractionBinsHF
int tower_diEta(int &iEta_1, int &iEta_2) const
std::vector< double > emFractionBinsBarrel
void SetJetClusterP4(double pt, double eta, double phi, double mass)
void produce(edm::Event &, const edm::EventSetup &) override
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
reco::Candidate::PolarLorentzVector seedTower
std::vector< double > tauCalibrationsBarrel
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< std::vector< std::vector< double > > > calibrationsHF
std::vector< std::vector< std::vector< std::vector< double > > > > tauPtCalibrationsHGCal
std::vector< double > jetPtBins
#define PI
Definition: QcdUeDQM.h:37
reco::Candidate::PolarLorentzVector l1egJetCluster
reco::Candidate::PolarLorentzVector jetCluster
void SetHcalJetClusterP4(double pt, double eta, double phi, double mass)
std::vector< double > tauL1egValuesBarrel
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
std::vector< std::vector< std::vector< std::vector< double > > > > tauPtCalibrationsBarrel
std::vector< edm::ParameterSet > tauL1egInfoHGCal
reco::Candidate::PolarLorentzVector GetP4() const
#define debug
Definition: HDRShower.cc:19
std::vector< double > tauPtBins
std::vector< double > absEtaBinsBarrel
double b
Definition: hdecay.h:120
L1CaloJetProducer(const edm::ParameterSet &)
std::vector< edm::ParameterSet > tauL1egInfoBarrel
std::vector< double > absEtaBinsHGCal
reco::Candidate::PolarLorentzVector ecalJetCluster
reco::Candidate::PolarLorentzVector p4
HLT enums.
double a
Definition: hdecay.h:121
std::vector< std::vector< std::vector< double > > > calibrationsHGCal
std::vector< double > jetCalibrationsHGCal
std::map< double, std::vector< double > > tauL1egInfoMapHGCal
edm::Handle< l1tp2::CaloTowerCollection > l1CaloTowerHandle
std::vector< double > tauCalibrationsHGCal
def move(src, dest)
Definition: eostools.py:511
std::vector< double > tauL1egValuesHGCal
#define LogDebug(id)
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:38