30 typedef ap_ufixed<14, 12, AP_RND, AP_WRAP>
pt_t;
32 const float ptLSB_ = 0.25;
33 const float phiLSB_ =
M_PI / 720;
36 typedef ap_fixed<pt_t::width + 1, pt_t::iwidth + 1, AP_RND, AP_SAT>
pxy_t;
37 typedef ap_fixed<2 * pt_t::width, 2 * pt_t::iwidth, AP_RND, AP_SAT>
pt2_t;
41 const phi_t hwPiOverTwo_ = round(
M_PI / (2 * phiLSB_));
43 typedef ap_ufixed<pt_t::width, 0>
inv_t;
46 const int dropBits_ = 2;
47 const int dropFactor_ = (1 << dropBits_);
48 const int invTableBits_ = 10;
49 const int invTableSize_ = (1 << invTableBits_);
59 maxCands_(
cfg.getParameter<
int>(
"maxCands")) {
60 produces<std::vector<l1t::EtSum>>();
67 std::vector<float>
pt;
68 std::vector<float>
phi;
71 const auto& l1PFCand = l1PFCandidates->at(
i);
72 pt.push_back(l1PFCand.pt());
73 phi.push_back(l1PFCand.phi());
82 std::unique_ptr<std::vector<l1t::EtSum>>
metCollection(
new std::vector<l1t::EtSum>(0));
88 std::vector<float>
phi,
99 Project(hw_pt, hw_phi, hw_px,
true);
100 Project(hw_pt, hw_phi, hw_py,
false);
102 hw_sumx = hw_sumx - hw_px;
103 hw_sumy = hw_sumy - hw_py;
107 hw_met =
sqrt(
int(hw_met));
109 phi_t hw_met_phi = 0;
112 metVector.SetPt(hw_met.to_double());
113 metVector.SetPhi(hw_met_phi.to_double() *
phiLSB_);
131 phiQ1 =
hwPi_ - phiQ1;
136 }
else if (phiQ1 < 0) {
141 typedef ap_ufixed<14, 12, AP_RND, AP_WRAP>
pt_t;
155 if (
px == 0 &&
py == 0) {
182 inv_t a_over_b =
a * inv_b;
185 printf(
" a, b = %f, %f; index, inv = %d, %f; ratio = %f \n",
190 a_over_b.to_double());
194 int atanTableBits_ = 7;
195 int atanTableSize_ = (1 << atanTableBits_);
196 index = round(a_over_b.to_double() * atanTableSize_);
200 printf(
" atan index, phi = %d, %f (%f rad) real atan(a/b)= %f \n",
204 atan(
a.to_double() /
b.to_double()));
211 if (px < 0 && py > 0)
213 if (
px > 0 &&
py < 0)
215 if (
px < 0 &&
py < 0)
219 printf(
" phi hw, float, real = %f, %f (%f rad from x,y = %f, %f) \n",
222 atan2(
py.to_double(),
px.to_double()),
~L1MetPfProducer() override
void CalcMetHLS(std::vector< float > pt, std::vector< float > phi, reco::Candidate::PolarLorentzVector &metVector) const
Sin< T >::type sin(const T &t)
void Project(pt_t pt, phi_t phi, pxy_t &pxy, bool isX, bool debug=false) const
ap_ufixed< 14, 12, AP_RND, AP_WRAP > pt_t
ap_fixed< 2 *pt_t::width, 2 *pt_t::iwidth, AP_RND, AP_SAT > pt2_t
void PhiFromXY(pxy_t px, pxy_t py, phi_t &phi, bool debug=false) const
Cos< T >::type cos(const T &t)
#define DEFINE_FWK_MODULE(type)
void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
ap_ufixed< pt_t::width, 0 > inv_t
edm::EDGetTokenT< vector< l1t::PFCandidate > > _l1PFToken
Log< level::Warning, false > LogWarning
L1MetPfProducer(const edm::ParameterSet &)
ap_fixed< pt_t::width+1, pt_t::iwidth+1, AP_RND, AP_SAT > pxy_t
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.