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