CMS 3D CMS Logo

L1GTScales.h
Go to the documentation of this file.
1 #ifndef L1Trigger_Phase2L1GT_L1GTScales_h
2 #define L1Trigger_Phase2L1GT_L1GTScales_h
3 
6 
7 #include <cmath>
8 
9 namespace l1t {
10  class L1GTScales {
11  public:
12  static constexpr int RELATIVE_ISOLATION_RESOLUTION = 18; // Resolution = 1/2^RELATIVE_ISOLATION_RESOLUTION
13  /* INV_MASS_SQR_OVER_2_DR_SQR_RESOLUTION originates from a simple analysis that yielded that the smallest
14  delta in hardware values of M^2/(2 dR^2) = 2.93326e-06 => log2(2.93326e-06) = -18.38 */
16 
17  L1GTScales(double pT_lsb,
18  double phi_lsb,
19  double eta_lsb,
20  double z0_lsb,
21  //double dD_lsb,
22  double isolationPT_lsb,
23  double beta_lsb,
24  double mass_lsb,
25  double seed_pT_lsb,
26  double seed_dZ_lsb,
27  double scalarSumPT_lsb,
28  double sum_pT_pv_lsb,
29  int pos_chg,
30  int neg_chg);
31 
33 
35 
36  int to_hw_pT_ceil(double value) const { return std::ceil(value / pT_lsb_); };
37  int to_hw_phi_ceil(double value) const { return std::ceil(value / phi_lsb_); };
38  int to_hw_eta_ceil(double value) const { return std::ceil(value / eta_lsb_); };
39  int to_hw_z0_ceil(double value) const { return std::ceil(value / z0_lsb_); };
40  // int to_hw_d0(double value) const { return std::ceil(value / d0_lsb_); };
43  }
44  int to_hw_isolationPT_ceil(double value) const { return std::ceil(value / isolationPT_lsb_); }
45  int to_hw_beta_ceil(double value) const { return std::ceil(value / beta_lsb_); };
46  int to_hw_mass_ceil(double value) const { return std::ceil(value / mass_lsb_); };
47  int to_hw_seed_pT_ceil(double value) const { return std::ceil(value / seed_pT_lsb_); };
48  int to_hw_seed_z0_ceil(double value) const { return std::ceil(value / seed_z0_lsb_); };
49  int to_hw_scalarSumPT_ceil(double value) const { return std::ceil(value / scalarSumPT_lsb_); };
50  int to_hw_sum_pT_pv_ceil(double value) const { return std::ceil(value / sum_pT_pv_lsb_); };
51 
52  int to_hw_dRSquared_ceil(double value) const { return std::ceil(value * value / (eta_lsb_ * eta_lsb_)); }
53 
54  double to_hw_InvMassSqrDiv2(double value) const { return value * value / (2 * pT_lsb_ * pT_lsb_); }
55  double to_hw_TransMassSqrDiv2(double value) const { return value * value / (2 * pT_lsb_ * pT_lsb_); }
56 
57  double to_hw_PtSquared(double value) const { return value * value / (pT_lsb_ * pT_lsb_); }
58  double to_hw_InvMassSqrOver2DR(double value) const {
60  (2 * pT_lsb_ * pT_lsb_);
61  }
62 
63  int to_hw_pT_floor(double value) const { return std::floor(value / pT_lsb_); };
64  int to_hw_phi_floor(double value) const { return std::floor(value / phi_lsb_); };
65  int to_hw_eta_floor(double value) const { return std::floor(value / eta_lsb_); };
66  int to_hw_z0_floor(double value) const { return std::floor(value / z0_lsb_); };
67  // int to_hw_d0(double value) const { return std::floor(value / d0_lsb_); };
70  }
71  int to_hw_isolationPT_floor(double value) const { return std::floor(value / isolationPT_lsb_); }
72  int to_hw_beta_floor(double value) const { return std::floor(value / beta_lsb_); };
73  int to_hw_mass_floor(double value) const { return std::floor(value / mass_lsb_); };
74  int to_hw_seed_pT_floor(double value) const { return std::floor(value / seed_pT_lsb_); };
75  int to_hw_seed_z0_floor(double value) const { return std::floor(value / seed_z0_lsb_); };
76  int to_hw_scalarSumPT_floor(double value) const { return std::floor(value / scalarSumPT_lsb_); };
77  int to_hw_sum_pT_pv_floor(double value) const { return std::floor(value / sum_pT_pv_lsb_); };
78 
79  int to_hw_dRSquared_floor(double value) const { return std::floor(value * value / (eta_lsb_ * eta_lsb_)); }
80 
81  double to_pT(int value) const { return value * pT_lsb_; };
82  double to_phi(int value) const { return value * phi_lsb_; };
83  double to_eta(int value) const { return value * eta_lsb_; };
84  double to_z0(int value) const { return value * z0_lsb_; };
85  double to_isolationPT(int value) const { return value * isolationPT_lsb_; }
86  double to_scalarSumPT(int value) const { return value * scalarSumPT_lsb_; };
87  int to_chg(int value) const { return value == pos_chg_ ? +1 : value == neg_chg_ ? -1 : 0; }
88 
89  double pT_lsb() const { return pT_lsb_; }
90  double phi_lsb() const { return phi_lsb_; }
91  double eta_lsb() const { return eta_lsb_; }
92  double z0_lsb() const { return z0_lsb_; }
93  double isolationPT_lsb() const { return isolationPT_lsb_; }
94  //const double dD_lsb_;
95  double beta_lsb() const { return beta_lsb_; }
96  double mass_lsb() const { return mass_lsb_; }
97  double seed_pT_lsb() const { return seed_pT_lsb_; }
98  double seed_z0_lsb() const { return seed_z0_lsb_; }
99  double scalarSumPT_lsb() const { return scalarSumPT_lsb_; }
100  double sum_pT_pv_lsb() const { return sum_pT_pv_lsb_; }
101  int pos_chg() const { return pos_chg_; }
102  int neg_chg() const { return neg_chg_; }
103 
104  private:
105  const double pT_lsb_;
106  const double phi_lsb_;
107  const double eta_lsb_;
108  const double z0_lsb_;
109  //const double dD_lsb_;
110  const double isolationPT_lsb_;
111  const double beta_lsb_;
112  const double mass_lsb_;
113  const double seed_pT_lsb_;
114  const double seed_z0_lsb_;
115  const double scalarSumPT_lsb_;
116  const double sum_pT_pv_lsb_;
117  const int pos_chg_;
118  const int neg_chg_;
119  };
120 } // namespace l1t
121 
122 #endif // L1Trigger_Phase2L1GT_L1GTScales_h
double to_hw_InvMassSqrOver2DR(double value) const
Definition: L1GTScales.h:58
double beta_lsb() const
Definition: L1GTScales.h:95
constexpr int32_t ceil(float num)
double eta_lsb() const
Definition: L1GTScales.h:91
double to_isolationPT(int value) const
Definition: L1GTScales.h:85
int to_hw_dRSquared_ceil(double value) const
Definition: L1GTScales.h:52
double pT_lsb() const
Definition: L1GTScales.h:89
int to_hw_mass_ceil(double value) const
Definition: L1GTScales.h:46
double to_hw_InvMassSqrDiv2(double value) const
Definition: L1GTScales.h:54
int to_hw_beta_ceil(double value) const
Definition: L1GTScales.h:45
double to_hw_TransMassSqrDiv2(double value) const
Definition: L1GTScales.h:55
int to_hw_z0_ceil(double value) const
Definition: L1GTScales.h:39
double phi_lsb() const
Definition: L1GTScales.h:90
int to_hw_sum_pT_pv_floor(double value) const
Definition: L1GTScales.h:77
int to_hw_scalarSumPT_ceil(double value) const
Definition: L1GTScales.h:49
const double seed_z0_lsb_
Definition: L1GTScales.h:114
int to_hw_phi_floor(double value) const
Definition: L1GTScales.h:64
int to_hw_seed_pT_ceil(double value) const
Definition: L1GTScales.h:47
int to_hw_isolationPT_floor(double value) const
Definition: L1GTScales.h:71
double to_phi(int value) const
Definition: L1GTScales.h:82
const double seed_pT_lsb_
Definition: L1GTScales.h:113
const double eta_lsb_
Definition: L1GTScales.h:107
int to_chg(int value) const
Definition: L1GTScales.h:87
static void fillPSetDescription(edm::ParameterSetDescription &)
Definition: L1GTScales.cc:53
delete x;
Definition: CaloConfig.h:22
int to_hw_beta_floor(double value) const
Definition: L1GTScales.h:72
double to_z0(int value) const
Definition: L1GTScales.h:84
int pos_chg() const
Definition: L1GTScales.h:101
double mass_lsb() const
Definition: L1GTScales.h:96
int to_hw_seed_pT_floor(double value) const
Definition: L1GTScales.h:74
int to_hw_sum_pT_pv_ceil(double value) const
Definition: L1GTScales.h:50
const double isolationPT_lsb_
Definition: L1GTScales.h:110
int to_hw_pT_floor(double value) const
Definition: L1GTScales.h:63
double to_eta(int value) const
Definition: L1GTScales.h:83
const int neg_chg_
Definition: L1GTScales.h:118
const double phi_lsb_
Definition: L1GTScales.h:106
int neg_chg() const
Definition: L1GTScales.h:102
const int pos_chg_
Definition: L1GTScales.h:117
int to_hw_relative_isolationPT_ceil(double value) const
Definition: L1GTScales.h:41
const double sum_pT_pv_lsb_
Definition: L1GTScales.h:116
int to_hw_dRSquared_floor(double value) const
Definition: L1GTScales.h:79
const double mass_lsb_
Definition: L1GTScales.h:112
int to_hw_seed_z0_ceil(double value) const
Definition: L1GTScales.h:48
Definition: value.py:1
static constexpr int INV_MASS_SQR_OVER_2_DR_SQR_RESOLUTION
Definition: L1GTScales.h:15
double scalarSumPT_lsb() const
Definition: L1GTScales.h:99
int to_hw_seed_z0_floor(double value) const
Definition: L1GTScales.h:75
const double scalarSumPT_lsb_
Definition: L1GTScales.h:115
const double beta_lsb_
Definition: L1GTScales.h:111
double sum_pT_pv_lsb() const
Definition: L1GTScales.h:100
double seed_pT_lsb() const
Definition: L1GTScales.h:97
int to_hw_mass_floor(double value) const
Definition: L1GTScales.h:73
double seed_z0_lsb() const
Definition: L1GTScales.h:98
int to_hw_eta_ceil(double value) const
Definition: L1GTScales.h:38
double to_pT(int value) const
Definition: L1GTScales.h:81
static constexpr int RELATIVE_ISOLATION_RESOLUTION
Definition: L1GTScales.h:12
const double pT_lsb_
Definition: L1GTScales.h:105
int to_hw_isolationPT_ceil(double value) const
Definition: L1GTScales.h:44
int to_hw_z0_floor(double value) const
Definition: L1GTScales.h:66
const double z0_lsb_
Definition: L1GTScales.h:108
int to_hw_pT_ceil(double value) const
Definition: L1GTScales.h:36
double z0_lsb() const
Definition: L1GTScales.h:92
int to_hw_eta_floor(double value) const
Definition: L1GTScales.h:65
double to_hw_PtSquared(double value) const
Definition: L1GTScales.h:57
int to_hw_relative_isolationPT_floor(double value) const
Definition: L1GTScales.h:68
double to_scalarSumPT(int value) const
Definition: L1GTScales.h:86
int to_hw_phi_ceil(double value) const
Definition: L1GTScales.h:37
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
int to_hw_scalarSumPT_floor(double value) const
Definition: L1GTScales.h:76
double isolationPT_lsb() const
Definition: L1GTScales.h:93
L1GTScales(double pT_lsb, double phi_lsb, double eta_lsb, double z0_lsb, double isolationPT_lsb, double beta_lsb, double mass_lsb, double seed_pT_lsb, double seed_dZ_lsb, double scalarSumPT_lsb, double sum_pT_pv_lsb, int pos_chg, int neg_chg)
Definition: L1GTScales.cc:8