CMS 3D CMS Logo

L1GTScales.cc
Go to the documentation of this file.
2 
3 #include <pybind11/pybind11.h>
4 
5 namespace py = pybind11;
6 
7 namespace l1t {
9  double phi_lsb,
10  double eta_lsb,
11  double z0_lsb,
12  //double d0_lsb,
13  double isolationPT_lsb,
14  double beta_lsb,
15  double mass_lsb,
16  double seed_pT_lsb,
17  double seed_z0_lsb,
18  double scalarSumPT_lsb,
19  double sum_pT_pv_lsb,
20  int pos_chg,
21  int neg_chg)
22  : pT_lsb_(pT_lsb),
23  phi_lsb_(phi_lsb),
24  eta_lsb_(eta_lsb),
25  z0_lsb_(z0_lsb),
26  //d0_lsb_(d0_lsb),
27  isolationPT_lsb_(isolationPT_lsb),
28  beta_lsb_(beta_lsb),
29  mass_lsb_(mass_lsb),
30  seed_pT_lsb_(seed_pT_lsb),
31  seed_z0_lsb_(seed_z0_lsb),
32  scalarSumPT_lsb_(scalarSumPT_lsb),
33  sum_pT_pv_lsb_(sum_pT_pv_lsb),
34  pos_chg_(pos_chg),
35  neg_chg_(neg_chg) {}
36 
38  : pT_lsb_(config.getParameter<double>("pT_lsb")),
39  phi_lsb_(config.getParameter<double>("phi_lsb")),
40  eta_lsb_(config.getParameter<double>("eta_lsb")),
41  z0_lsb_(config.getParameter<double>("z0_lsb")),
42  //d0_lsb_(config.getParameter<double>("d0_lsb")),
43  isolationPT_lsb_(config.getParameter<double>("isolationPT_lsb")),
44  beta_lsb_(config.getParameter<double>("beta_lsb")),
45  mass_lsb_(config.getParameter<double>("mass_lsb")),
46  seed_pT_lsb_(config.getParameter<double>("seed_pT_lsb")),
47  seed_z0_lsb_(config.getParameter<double>("seed_z0_lsb")),
48  scalarSumPT_lsb_(config.getParameter<double>("scalarSumPT_lsb")),
49  sum_pT_pv_lsb_(config.getParameter<double>("sum_pT_pv_lsb")),
50  pos_chg_(config.getParameter<int>("pos_chg")),
51  neg_chg_(config.getParameter<int>("neg_chg")) {}
52 
54  desc.add<double>("pT_lsb");
55  desc.add<double>("phi_lsb");
56  desc.add<double>("eta_lsb");
57  desc.add<double>("z0_lsb");
58  //desc.add<double>("d0_lsb");
59  desc.add<double>("isolationPT_lsb");
60  desc.add<double>("beta_lsb");
61  desc.add<double>("mass_lsb");
62  desc.add<double>("seed_pT_lsb");
63  desc.add<double>("seed_z0_lsb");
64  desc.add<double>("scalarSumPT_lsb");
65  desc.add<double>("sum_pT_pv_lsb");
66  desc.add<int>("pos_chg");
67  desc.add<int>("neg_chg");
68  }
69 
70  PYBIND11_MODULE(libL1TriggerPhase2L1GT, m) {
71  py::class_<L1GTScales>(m, "L1GTScales")
72  .def(py::init<double,
73  double,
74  double,
75  double,
76  /*double, */
77  double,
78  double,
79  double,
80  double,
81  double,
82  double,
83  double,
84  int,
85  int>())
86  .def("to_hw_pT_ceil", &L1GTScales::to_hw_pT_ceil)
87  .def("to_hw_phi_ceil", &L1GTScales::to_hw_phi_ceil)
88  .def("to_hw_eta_ceil", &L1GTScales::to_hw_eta_ceil)
89  .def("to_hw_z0_ceil", &L1GTScales::to_hw_z0_ceil)
90  .def("to_hw_relative_isolationPT_ceil", &L1GTScales::to_hw_relative_isolationPT_ceil)
91  .def("to_hw_beta_ceil", &L1GTScales::to_hw_beta_ceil)
92  .def("to_hw_mass_ceil", &L1GTScales::to_hw_mass_ceil)
93  .def("to_hw_seed_pT_ceil", &L1GTScales::to_hw_seed_pT_ceil)
94  .def("to_hw_seed_z0_ceil", &L1GTScales::to_hw_seed_z0_ceil)
95  .def("to_hw_scalarSumPT_ceil", &L1GTScales::to_hw_scalarSumPT_ceil)
96  .def("to_hw_sum_pT_pv_ceil", &L1GTScales::to_hw_sum_pT_pv_ceil)
97  .def("to_hw_dRSquared_ceil", &L1GTScales::to_hw_dRSquared_ceil)
98  .def("to_hw_pT_floor", &L1GTScales::to_hw_pT_floor)
99  .def("to_hw_phi_floor", &L1GTScales::to_hw_phi_floor)
100  .def("to_hw_eta_floor", &L1GTScales::to_hw_eta_floor)
101  .def("to_hw_z0_floor", &L1GTScales::to_hw_z0_floor)
102  .def("to_hw_relative_isolationPT_floor", &L1GTScales::to_hw_relative_isolationPT_floor)
103  .def("to_hw_beta_floor", &L1GTScales::to_hw_beta_floor)
104  .def("to_hw_mass_floor", &L1GTScales::to_hw_mass_floor)
105  .def("to_hw_seed_pT_floor", &L1GTScales::to_hw_seed_pT_floor)
106  .def("to_hw_seed_z0_floor", &L1GTScales::to_hw_seed_z0_floor)
107  .def("to_hw_scalarSumPT_floor", &L1GTScales::to_hw_scalarSumPT_floor)
108  .def("to_hw_sum_pT_pv_floor", &L1GTScales::to_hw_sum_pT_pv_floor)
109  .def("to_hw_dRSquared_floor", &L1GTScales::to_hw_dRSquared_floor)
110  .def("to_hw_InvMassSqrDiv2", &L1GTScales::to_hw_InvMassSqrDiv2)
111  .def("to_hw_TransMassSqrDiv2", &L1GTScales::to_hw_TransMassSqrDiv2)
112  .def("to_hw_PtSquared", &L1GTScales::to_hw_PtSquared)
113  .def("to_hw_InvMassSqrOver2DR", &L1GTScales::to_hw_InvMassSqrOver2DR)
114  .def("neg_chg", &L1GTScales::neg_chg)
115  .def("pos_chg", &L1GTScales::pos_chg);
116  }
117 } // namespace l1t
double to_hw_InvMassSqrOver2DR(double value) const
Definition: L1GTScales.h:58
int def(FILE *, FILE *, int)
int to_hw_dRSquared_ceil(double value) const
Definition: L1GTScales.h:52
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
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
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 init
Definition: HydjetWrapper.h:66
static void fillPSetDescription(edm::ParameterSetDescription &)
Definition: L1GTScales.cc:53
delete x;
Definition: CaloConfig.h:22
Definition: config.py:1
int to_hw_beta_floor(double value) const
Definition: L1GTScales.h:72
int pos_chg() const
Definition: L1GTScales.h:101
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
int to_hw_pT_floor(double value) const
Definition: L1GTScales.h:63
int neg_chg() const
Definition: L1GTScales.h:102
int to_hw_relative_isolationPT_ceil(double value) const
Definition: L1GTScales.h:41
int to_hw_dRSquared_floor(double value) const
Definition: L1GTScales.h:79
int to_hw_seed_z0_ceil(double value) const
Definition: L1GTScales.h:48
int to_hw_seed_z0_floor(double value) const
Definition: L1GTScales.h:75
int to_hw_mass_floor(double value) const
Definition: L1GTScales.h:73
int to_hw_eta_ceil(double value) const
Definition: L1GTScales.h:38
int to_hw_z0_floor(double value) const
Definition: L1GTScales.h:66
PYBIND11_MODULE(libL1TriggerPhase2L1GT, m)
Definition: L1GTScales.cc:70
int to_hw_pT_ceil(double value) const
Definition: L1GTScales.h:36
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
int to_hw_phi_ceil(double value) const
Definition: L1GTScales.h:37
int to_hw_scalarSumPT_floor(double value) const
Definition: L1GTScales.h:76
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