CMS 3D CMS Logo

L1MetPfProducer.cc
Go to the documentation of this file.
1 #include <vector>
2 #include <ap_int.h>
3 #include <ap_fixed.h>
4 #include <TVector2.h>
5 
11 
15 
16 using namespace l1t;
17 
19 public:
20  explicit L1MetPfProducer(const edm::ParameterSet&);
21  ~L1MetPfProducer() override;
22 
23 private:
24  void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
26 
27  int maxCands_ = 128;
28 
29  // quantization controllers
30  typedef ap_ufixed<14, 12, AP_RND, AP_WRAP> pt_t; // LSB is 0.25 and max is 4 TeV
31  typedef ap_int<12> phi_t; // LSB is pi/720 ~ 0.0044 and max is +/-8.9
32  const float ptLSB_ = 0.25; // GeV
33  const float phiLSB_ = M_PI / 720; // rad
34 
35  // derived, helper types
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;
38  // derived, helper constants
39  const float maxPt_ = ((1 << pt_t::width) - 1) * ptLSB_;
40  const phi_t hwPi_ = round(M_PI / phiLSB_);
41  const phi_t hwPiOverTwo_ = round(M_PI / (2 * phiLSB_));
42 
43  typedef ap_ufixed<pt_t::width, 0> inv_t; // can't easily use the MAXPT/pt trick with ap_fixed
44 
45  // to make configurable...
46  const int dropBits_ = 2;
47  const int dropFactor_ = (1 << dropBits_);
48  const int invTableBits_ = 10;
49  const int invTableSize_ = (1 << invTableBits_);
50 
51  void Project(pt_t pt, phi_t phi, pxy_t& pxy, bool isX, bool debug = false) const;
52  void PhiFromXY(pxy_t px, pxy_t py, phi_t& phi, bool debug = false) const;
53 
54  void CalcMetHLS(std::vector<float> pt, std::vector<float> phi, reco::Candidate::PolarLorentzVector& metVector) const;
55 };
56 
58  : _l1PFToken(consumes<std::vector<l1t::PFCandidate>>(cfg.getParameter<edm::InputTag>("L1PFObjects"))),
59  maxCands_(cfg.getParameter<int>("maxCands")) {
60  produces<std::vector<l1t::EtSum>>();
61 }
62 
65  iEvent.getByToken(_l1PFToken, l1PFCandidates);
66 
67  std::vector<float> pt;
68  std::vector<float> phi;
69 
70  for (int i = 0; i < int(l1PFCandidates->size()) && (i < maxCands_ || maxCands_ < 0); i++) {
71  const auto& l1PFCand = l1PFCandidates->at(i);
72  pt.push_back(l1PFCand.pt());
73  phi.push_back(l1PFCand.phi());
74  }
75 
77 
78  CalcMetHLS(pt, phi, metVector);
79 
80  l1t::EtSum theMET(metVector, l1t::EtSum::EtSumType::kTotalHt, 0, 0, 0, 0);
81 
82  std::unique_ptr<std::vector<l1t::EtSum>> metCollection(new std::vector<l1t::EtSum>(0));
83  metCollection->push_back(theMET);
85 }
86 
87 void L1MetPfProducer::CalcMetHLS(std::vector<float> pt,
88  std::vector<float> phi,
89  reco::Candidate::PolarLorentzVector& metVector) const {
90  pxy_t hw_px = 0;
91  pxy_t hw_py = 0;
92  pxy_t hw_sumx = 0;
93  pxy_t hw_sumy = 0;
94 
95  for (uint i = 0; i < pt.size(); i++) {
96  pt_t hw_pt = min(pt[i], maxPt_);
97  phi_t hw_phi = float(TVector2::Phi_mpi_pi(phi[i]) / phiLSB_);
98 
99  Project(hw_pt, hw_phi, hw_px, true);
100  Project(hw_pt, hw_phi, hw_py, false);
101 
102  hw_sumx = hw_sumx - hw_px;
103  hw_sumy = hw_sumy - hw_py;
104  }
105 
106  pt2_t hw_met = pt2_t(hw_sumx) * pt2_t(hw_sumx) + pt2_t(hw_sumy) * pt2_t(hw_sumy);
107  hw_met = sqrt(int(hw_met)); // stand-in for HLS::sqrt
108 
109  phi_t hw_met_phi = 0;
110  PhiFromXY(hw_sumx, hw_sumy, hw_met_phi);
111 
112  metVector.SetPt(hw_met.to_double());
113  metVector.SetPhi(hw_met_phi.to_double() * phiLSB_);
114  metVector.SetEta(0);
115 }
116 
117 void L1MetPfProducer::Project(pt_t pt, phi_t phi, pxy_t& pxy, bool isX, bool debug) const {
118  /*
119  Convert pt and phi to px (py)
120  1) Map phi to the first quadrant to reduce LUT size
121  2) Lookup sin(phiQ1), where the result is in [0,maxPt]
122  which is used to encode [0,1].
123  3) Multiply pt by sin(phiQ1) to get px. Result will be px*maxPt, but
124  wrapping multiplication is 'mod maxPt' so the correct value is returned.
125  4) Check px=-|px|.
126  */
127 
128  // set phi to first quadrant
129  phi_t phiQ1 = (phi > 0) ? phi : phi_t(-phi); // Q1/Q4
130  if (phiQ1 >= hwPiOverTwo_)
131  phiQ1 = hwPi_ - phiQ1;
132 
133  if (phiQ1 > hwPiOverTwo_) {
134  edm::LogWarning("L1MetPfProducer") << "unexpected phi (high)";
135  phiQ1 = hwPiOverTwo_;
136  } else if (phiQ1 < 0) {
137  edm::LogWarning("L1MetPfProducer") << "unexpected phi (low)";
138  phiQ1 = 0;
139  }
140  if (isX) {
141  typedef ap_ufixed<14, 12, AP_RND, AP_WRAP> pt_t; // LSB is 0.25 and max is 4 TeV
142  ap_ufixed<pt_t::width, 0> cosPhi = cos(phiQ1.to_double() / hwPiOverTwo_.to_double() * M_PI / 2);
143  pxy = pt * cosPhi;
144  if (phi > hwPiOverTwo_ || phi < -hwPiOverTwo_)
145  pxy = -pxy;
146  } else {
147  ap_ufixed<pt_t::width, 0> sinPhi = sin(phiQ1.to_double() / hwPiOverTwo_.to_double() * M_PI / 2);
148  pxy = pt * sinPhi;
149  if (phi < 0)
150  pxy = -pxy;
151  }
152 }
153 
155  if (px == 0 && py == 0) {
156  phi = 0;
157  return;
158  }
159  if (px == 0) {
160  phi = py > 0 ? hwPiOverTwo_ : phi_t(-hwPiOverTwo_);
161  return;
162  }
163  if (py == 0) {
164  phi = px > 0 ? phi_t(0) : phi_t(-hwPi_);
165  return;
166  }
167 
168  // get q1 coordinates
169  pt_t x = px > 0 ? pt_t(px) : pt_t(-px); //px>=0 ? px : -px;
170  pt_t y = py > 0 ? pt_t(py) : pt_t(-py); //px>=0 ? px : -px;
171  // transform so a<b
172  pt_t a = x < y ? x : y;
173  pt_t b = x < y ? y : x;
174 
175  if (b.to_double() > maxPt_ / dropFactor_)
176  b = maxPt_ / dropFactor_;
177  // map [0,max/4) to inv table size
178  int index = round((b.to_double() / (maxPt_ / dropFactor_)) * invTableSize_);
179  float bcheck = (float(index) / invTableSize_) * (maxPt_ / dropFactor_);
180  inv_t inv_b = 1. / ((float(index) / invTableSize_) * (maxPt_ / dropFactor_));
181 
182  inv_t a_over_b = a * inv_b;
183 
184  if (debug) {
185  printf(" a, b = %f, %f; index, inv = %d, %f; ratio = %f \n",
186  a.to_double(),
187  b.to_double(),
188  index,
189  inv_b.to_double(),
190  a_over_b.to_double());
191  printf(" bcheck, 1/bc = %f, %f -- %d %f %d \n", bcheck, 1. / bcheck, invTableSize_, maxPt_, dropFactor_);
192  }
193 
194  int atanTableBits_ = 7;
195  int atanTableSize_ = (1 << atanTableBits_);
196  index = round(a_over_b.to_double() * atanTableSize_);
197  phi = atan(float(index) / atanTableSize_) / phiLSB_;
198 
199  if (debug) {
200  printf(" atan index, phi = %d, %f (%f rad) real atan(a/b)= %f \n",
201  index,
202  phi.to_double(),
203  phi.to_double() * (M_PI / hwPi_.to_double()),
204  atan(a.to_double() / b.to_double()));
205  }
206 
207  // rotate from (0,pi/4) to full quad1
208  if (y > x)
209  phi = hwPiOverTwo_ - phi; //phi = pi/2 - phi
210  // other quadrants
211  if (px < 0 && py > 0)
212  phi = hwPi_ - phi; // Q2 phi = pi - phi
213  if (px > 0 && py < 0)
214  phi = -phi; // Q4 phi = -phi
215  if (px < 0 && py < 0)
216  phi = -(hwPi_ - phi); // Q3 composition of both
217 
218  if (debug) {
219  printf(" phi hw, float, real = %f, %f (%f rad from x,y = %f, %f) \n",
220  phi.to_double(),
221  phi.to_double() * (M_PI / hwPi_.to_double()),
222  atan2(py.to_double(), px.to_double()),
223  px.to_double(),
224  py.to_double());
225  }
226 }
227 
229 
~L1MetPfProducer() override
const int dropFactor_
void CalcMetHLS(std::vector< float > pt, std::vector< float > phi, reco::Candidate::PolarLorentzVector &metVector) const
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
delete x;
Definition: CaloConfig.h:22
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
const float maxPt_
int iEvent
Definition: GenABIO.cc:224
void PhiFromXY(pxy_t px, pxy_t py, phi_t &phi, bool debug=false) const
const phi_t hwPi_
T sqrt(T t)
Definition: SSEVec.h:19
const phi_t hwPiOverTwo_
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
ap_int< 12 > phi_t
#define M_PI
ap_ufixed< pt_t::width, 0 > inv_t
#define debug
Definition: HDRShower.cc:19
double b
Definition: hdecay.h:120
HLT enums.
double a
Definition: hdecay.h:121
edm::EDGetTokenT< vector< l1t::PFCandidate > > _l1PFToken
Log< level::Warning, false > LogWarning
L1MetPfProducer(const edm::ParameterSet &)
const float phiLSB_
const int invTableSize_
ap_fixed< pt_t::width+1, pt_t::iwidth+1, AP_RND, AP_SAT > pxy_t
def move(src, dest)
Definition: eostools.py:511
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:38