CMS 3D CMS Logo

L1TCaloLayer1FetchLUTs.cc
Go to the documentation of this file.
1 // This function fetches Layer1 ECAL and HCAL LUTs from CMSSW configuration
2 // It is provided as a global helper function outside of class structure
3 // so that it can be shared by L1CaloLayer1 and L1CaloLayer1Spy
4 
5 #include <vector>
6 
10 
14 
20 
23 
24 #include "L1TCaloLayer1FetchLUTs.hh"
25 #include "UCTLogging.hh"
26 
27 using namespace l1tcalo;
28 
30  const L1TCaloLayer1FetchLUTsTokens &iTokens,
31  const edm::EventSetup &iSetup,
32  std::vector<std::array<std::array<std::array<uint32_t, nEtBins>, nCalSideBins>, nCalEtaBins> > &eLUT,
33  std::vector<std::array<std::array<std::array<uint32_t, nEtBins>, nCalSideBins>, nCalEtaBins> > &hLUT,
34  std::vector<std::array<std::array<uint32_t, nEtBins>, nHfEtaBins> > &hfLUT,
35  std::vector<unsigned long long int> &hcalFBLUT,
36  std::vector<unsigned int> &ePhiMap,
37  std::vector<unsigned int> &hPhiMap,
38  std::vector<unsigned int> &hfPhiMap,
39  bool useLSB,
40  bool useCalib,
41  bool useECALLUT,
42  bool useHCALLUT,
43  bool useHFLUT,
44  bool useHCALFBLUT,
45  int fwVersion) {
46  int hfValid = 1;
47  const HcalTrigTowerGeometry &pG = iSetup.getData(iTokens.geom_);
48  if (!pG.use1x1()) {
49  edm::LogError("L1TCaloLayer1FetchLUTs")
50  << "Using Stage2-Layer1 but HCAL Geometry has use1x1 = 0! HF will be suppressed. Check Global Tag, etc.";
51  hfValid = 0;
52  }
53 
54  // CaloParams contains all persisted parameters for Layer 1
55  edm::ESHandle<l1t::CaloParams> paramsHandle = iSetup.getHandle(iTokens.params_);
56  if (not paramsHandle.isValid()) {
57  edm::LogError("L1TCaloLayer1FetchLUTs") << "Missing CaloParams object! Check Global Tag, etc.";
58  return false;
59  }
60  l1t::CaloParamsHelper caloParams(*paramsHandle.product());
61 
62  // Calo Trigger Layer1 output LSB Real ET value
63  double caloLSB = caloParams.towerLsbSum();
64  if (caloLSB != 0.5) {
65  // Lots of things expect this, better give fair warning if not
66  edm::LogError("L1TCaloLayer1FetchLUTs") << "caloLSB (caloParams.towerLsbSum()) != 0.5, actually = " << caloLSB;
67  }
68 
69  // ECal/HCal scale factors will be a x*y*28 array:
70  // ieta = 28 eta scale factors (1 .. 28)
71  // etBin = size of Real ET Bins vector
72  // phiBin = max(Real Phi Bins vector)
73  // So, index = phiBin*etBin*28+etBin*28+ieta
74  auto ecalScaleETBins = caloParams.layer1ECalScaleETBins();
75  auto ecalScalePhiBins = caloParams.layer1ECalScalePhiBins();
76  if (ecalScalePhiBins.empty()) {
77  // Backwards-compatibility (no phi binning)
78  ecalScalePhiBins.resize(36, 0);
79  } else if (ecalScalePhiBins.size() % 36 != 0) {
80  edm::LogError("L1TCaloLayer1FetchLUTs") << "caloParams.layer1ECalScaleETBins().size() is not multiple of 36 !!";
81  return false;
82  }
83  size_t numEcalPhiBins = (*std::max_element(ecalScalePhiBins.begin(), ecalScalePhiBins.end())) + 1;
84  auto ecalSF = caloParams.layer1ECalScaleFactors();
85  if (ecalSF.size() != ecalScaleETBins.size() * numEcalPhiBins * 28) {
86  edm::LogError("L1TCaloLayer1FetchLUTs") << "caloParams.layer1ECalScaleFactors().size() != "
87  "caloParams.layer1ECalScaleETBins().size()*numEcalPhiBins*28 !!";
88  return false;
89  }
90  auto hcalScaleETBins = caloParams.layer1HCalScaleETBins();
91  auto hcalScalePhiBins = caloParams.layer1HCalScalePhiBins();
92  if (hcalScalePhiBins.empty()) {
93  hcalScalePhiBins.resize(36, 0);
94  } else if (hcalScalePhiBins.size() % 36 != 0) {
95  edm::LogError("L1TCaloLayer1FetchLUTs") << "caloParams.layer1HCalScaleETBins().size() is not multiple of 36 !!";
96  return false;
97  }
98  size_t numHcalPhiBins = (*std::max_element(hcalScalePhiBins.begin(), hcalScalePhiBins.end())) + 1;
99  auto hcalSF = caloParams.layer1HCalScaleFactors();
100  if (hcalSF.size() != hcalScaleETBins.size() * numHcalPhiBins * 28) {
101  edm::LogError("L1TCaloLayer1FetchLUTs") << "caloParams.layer1HCalScaleFactors().size() != "
102  "caloParams.layer1HCalScaleETBins().size()*numHcalPhiBins*28 !!";
103  return false;
104  }
105 
106  // HF 1x1 scale factors will be a x*y*12 array:
107  // ieta = 12 eta scale factors (30 .. 41)
108  // etBin = size of Real ET Bins vector
109  // phiBin = max(Real Phi Bins vector)
110  // So, index = phiBin*etBin*12+etBin*12+ieta
111  auto hfScaleETBins = caloParams.layer1HFScaleETBins();
112  auto hfScalePhiBins = caloParams.layer1HFScalePhiBins();
113  if (hfScalePhiBins.empty()) {
114  hfScalePhiBins.resize(36, 0);
115  } else if (hfScalePhiBins.size() % 36 != 0) {
116  edm::LogError("L1TCaloLayer1FetchLUTs") << "caloParams.layer1HFScaleETBins().size() is not multiple of 36 !!";
117  return false;
118  }
119  size_t numHFPhiBins = (*std::max_element(hfScalePhiBins.begin(), hfScalePhiBins.end())) + 1;
120  auto hfSF = caloParams.layer1HFScaleFactors();
121  if (hfSF.size() != hfScaleETBins.size() * numHFPhiBins * 12) {
122  edm::LogError("L1TCaloLayer1FetchLUTs")
123  << "caloParams.layer1HFScaleFactors().size() != caloParams.layer1HFScaleETBins().size()*numHFPhiBins*12 !!";
124  return false;
125  }
126 
127  // Sanity check scale factors exist
128  if (useCalib && (ecalSF.empty() || hcalSF.empty() || hfSF.empty())) {
129  edm::LogError("L1TCaloLayer1FetchLUTs") << "Layer 1 calibrations requested (useCalib = True) but there are missing "
130  "scale factors in CaloParams! Please check conditions setup.";
131  return false;
132  }
133 
134  // HCAL FB LUT will be a 1*28 array:
135  // ieta = 28 eta scale factors (1 .. 28)
136  // So, index = ieta
137  auto fbLUTUpper = caloParams.layer1HCalFBLUTUpper();
138  auto fbLUTLower = caloParams.layer1HCalFBLUTLower();
139  // Only check for HCAL FB LUT if useHCALFBLUT = true
140  if (useHCALFBLUT) {
141  if (fbLUTUpper.size() != nCalEtaBins) {
142  edm::LogError("L1TCaloLayer1FetchLUTs")
143  << "caloParams.layer1HCalFBLUTUpper().size() " << fbLUTUpper.size() << " != " << nCalEtaBins << " !!";
144  return false;
145  }
146  if (fbLUTLower.size() != nCalEtaBins) {
147  edm::LogError("L1TCaloLayer1FetchLUTs")
148  << "caloParams.layer1HCalFBLUTLower().size() " << fbLUTLower.size() << " != " << nCalEtaBins << " !!";
149  return false;
150  }
151  }
152 
153  // get energy scale to convert input from ECAL - this should be linear with LSB = 0.5 GeV
154  const double ecalLSB = 0.5;
155 
156  // get energy scale to convert input from HCAL
157  edm::ESHandle<CaloTPGTranscoder> decoder = iSetup.getHandle(iTokens.decoder_);
158  if (not decoder.isValid()) {
159  edm::LogError("L1TCaloLayer1FetchLUTs") << "Missing CaloTPGTranscoder object! Check Global Tag, etc.";
160  return false;
161  }
162 
163  // TP compression scale is always phi symmetric
164  // We default to 3 since HF has no ieta=41 iphi=1,2
165  auto decodeHcalEt = [&decoder](int iEta, uint32_t compressedEt, uint32_t iPhi = 3) -> double {
166  HcalTriggerPrimitiveSample sample(compressedEt);
167  HcalTrigTowerDetId id(iEta, iPhi);
168  if (std::abs(iEta) >= 30) {
169  id.setVersion(1);
170  }
171  return decoder->hcaletValue(id, sample);
172  };
173 
174  // Make ECal LUT
175  for (uint32_t phiBin = 0; phiBin < numEcalPhiBins; phiBin++) {
176  std::array<std::array<std::array<uint32_t, nEtBins>, nCalSideBins>, nCalEtaBins> phiLUT;
177  eLUT.push_back(phiLUT);
178  for (uint32_t etaBin = 0; etaBin < nCalEtaBins; etaBin++) {
179  for (uint32_t fb = 0; fb < nCalSideBins; fb++) {
180  for (uint32_t ecalInput = 0; ecalInput <= 0xFF; ecalInput++) {
181  uint32_t value = ecalInput;
182  if (useECALLUT) {
183  double linearizedECalInput = ecalInput * ecalLSB; // in GeV
184 
185  uint32_t etBin = 0;
186  for (; etBin < ecalScaleETBins.size(); etBin++) {
187  if (linearizedECalInput < ecalScaleETBins[etBin])
188  break;
189  }
190  if (etBin >= ecalScaleETBins.size())
191  etBin = ecalScaleETBins.size() - 1;
192 
193  double calibratedECalInput = linearizedECalInput;
194  if (useCalib)
195  calibratedECalInput *= ecalSF.at(phiBin * ecalScaleETBins.size() * 28 + etBin * 28 + etaBin);
196  if (useLSB)
197  calibratedECalInput /= caloLSB;
198 
199  value = calibratedECalInput;
200  if (fwVersion > 2) {
201  // Saturate if either decompressed value is over 127.5 GeV or input saturated
202  // (meaningless for ecal, since ecalLSB == caloLSB)
203  if (value > 0xFF || ecalInput == 0xFF) {
204  value = 0xFF;
205  }
206  } else {
207  if (value > 0xFF) {
208  value = 0xFF;
209  }
210  }
211  }
212  if (value == 0) {
213  value = (1 << 11);
214  } else {
215  uint32_t et_log2 = ((uint32_t)log2(value)) & 0x7;
216  value |= (et_log2 << 12);
217  }
218  value |= (fb << 10);
219  eLUT[phiBin][etaBin][fb][ecalInput] = value;
220  }
221  }
222  }
223  }
224 
225  // Make HCal LUT
226  for (uint32_t phiBin = 0; phiBin < numHcalPhiBins; phiBin++) {
227  std::array<std::array<std::array<uint32_t, nEtBins>, nCalSideBins>, nCalEtaBins> phiLUT;
228  hLUT.push_back(phiLUT);
229  for (uint32_t etaBin = 0; etaBin < nCalEtaBins; etaBin++) {
230  int caloEta = etaBin + 1;
231  int iPhi = 3;
232  auto pos = std::find(hcalScalePhiBins.begin(), hcalScalePhiBins.end(), phiBin);
233  if (pos != hcalScalePhiBins.end()) {
234  // grab an iPhi bin
235  auto index = std::distance(hcalScalePhiBins.begin(), pos);
236  if (index < 18) {
237  caloEta *= -1;
238  iPhi = index * 4 + 1;
239  } else {
240  iPhi = (index - 18) * 4 + 1;
241  }
242  }
243  for (uint32_t fb = 0; fb < nCalSideBins; fb++) {
244  for (uint32_t hcalInput = 0; hcalInput <= 0xFF; hcalInput++) {
245  uint32_t value = hcalInput;
246  if (useHCALLUT) {
247  // hcaletValue defined in L137 of CalibCalorimetry/CaloTPG/src/CaloTPGTranscoderULUT.cc
248  double linearizedHcalInput = decodeHcalEt(caloEta, hcalInput, iPhi); // in GeV
249 
250  uint32_t etBin = 0;
251  for (; etBin < hcalScaleETBins.size(); etBin++) {
252  if (linearizedHcalInput < hcalScaleETBins[etBin])
253  break;
254  }
255  if (etBin >= hcalScaleETBins.size())
256  etBin = hcalScaleETBins.size() - 1;
257 
258  double calibratedHcalInput = linearizedHcalInput;
259  if (useCalib)
260  calibratedHcalInput *= hcalSF.at(phiBin * hcalScaleETBins.size() * 28 + etBin * 28 + etaBin);
261  if (useLSB)
262  calibratedHcalInput /= caloLSB;
263 
264  value = calibratedHcalInput;
265  if (fwVersion > 2) {
266  // Saturate if either decompressed value is over 127.5 GeV or input saturated
267  if (value > 0xFF || hcalInput == 0xFF) {
268  value = 0xFF;
269  }
270  } else {
271  if (value > 0xFF) {
272  value = 0xFF;
273  }
274  }
275  }
276  if (value == 0) {
277  value = (1 << 11);
278  } else {
279  uint32_t et_log2 = ((uint32_t)log2(value)) & 0x7;
280  value |= (et_log2 << 12);
281  }
282  value |= (fb << 10);
283  hLUT[phiBin][etaBin][fb][hcalInput] = value;
284  }
285  }
286  }
287  }
288 
289  // Make HF LUT
290  for (uint32_t phiBin = 0; phiBin < numHFPhiBins; phiBin++) {
291  std::array<std::array<uint32_t, nEtBins>, nHfEtaBins> phiLUT;
292  hfLUT.push_back(phiLUT);
293  for (uint32_t etaBin = 0; etaBin < nHfEtaBins; etaBin++) {
294  int caloEta = etaBin + 30;
295  int iPhi = 3;
296  auto pos = std::find(hfScalePhiBins.begin(), hfScalePhiBins.end(), phiBin);
297  if (pos != hfScalePhiBins.end()) {
298  auto index = std::distance(hfScalePhiBins.begin(), pos);
299  if (index < 18) {
300  caloEta *= -1;
301  iPhi = index * 4 - 1;
302  } else {
303  iPhi = (index - 18) * 4 - 1;
304  }
305  if (iPhi < 0)
306  iPhi = 71;
307  }
308  for (uint32_t etCode = 0; etCode < nEtBins; etCode++) {
309  uint32_t value = etCode;
310  if (useHFLUT) {
311  double linearizedHFInput = 0;
312  if (hfValid) {
313  linearizedHFInput = decodeHcalEt(caloEta, value, iPhi); // in GeV
314  }
315 
316  uint32_t etBin = 0;
317  for (; etBin < hfScaleETBins.size(); etBin++) {
318  if (linearizedHFInput < hfScaleETBins[etBin])
319  break;
320  }
321  if (etBin >= hfScaleETBins.size())
322  etBin = hfScaleETBins.size() - 1;
323 
324  double calibratedHFInput = linearizedHFInput;
325  if (useCalib)
326  calibratedHFInput *= hfSF.at(phiBin * hfScalePhiBins.size() * 12 + etBin * 12 + etaBin);
327  if (useLSB)
328  calibratedHFInput /= caloLSB;
329 
330  if (fwVersion > 2) {
331  uint32_t absCaloEta = std::abs(caloEta);
332  if (absCaloEta > 29 && absCaloEta < 40) {
333  // Divide by two (since two duplicate towers are sent)
334  calibratedHFInput *= 0.5;
335  } else if (absCaloEta == 40 || absCaloEta == 41) {
336  // Divide by four
337  calibratedHFInput *= 0.25;
338  }
339  value = calibratedHFInput;
340  // Saturate if either decompressed value is over 127.5 GeV or input saturated
341  if (value >= 0xFF || etCode == 0xFF) {
342  value = 0x1FD;
343  }
344  } else {
345  value = calibratedHFInput;
346  if (value > 0xFF) {
347  value = 0xFF;
348  }
349  }
350  }
351  hfLUT[phiBin][etaBin][etCode] = value;
352  }
353  }
354  }
355 
356  // Make HCal FB LUT
357  for (uint32_t etaBin = 0; etaBin < nCalEtaBins; etaBin++) {
358  uint64_t value = 0xFFFFFFFFFFFFFFFF;
359  if (useHCALFBLUT) {
360  value = (((uint64_t)fbLUTUpper.at(etaBin)) << 32) | fbLUTLower.at(etaBin);
361  }
362  hcalFBLUT.push_back(value);
363  }
364 
365  // plus/minus, 18 CTP7, 4 iPhi each
366  for (uint32_t isPos = 0; isPos < 2; isPos++) {
367  for (uint32_t iPhi = 1; iPhi <= 72; iPhi++) {
368  uint32_t card = floor((iPhi + 1) / 4);
369  if (card > 17)
370  card -= 18;
371  ePhiMap[isPos * 72 + iPhi - 1] = ecalScalePhiBins[isPos * 18 + card];
372  hPhiMap[isPos * 72 + iPhi - 1] = hcalScalePhiBins[isPos * 18 + card];
373  hfPhiMap[isPos * 72 + iPhi - 1] = hfScalePhiBins[isPos * 18 + card];
374  }
375  }
376 
377  return true;
378 }
379 /* vim: set ts=8 sw=2 tw=0 et :*/
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
Log< level::Error, false > LogError
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
T const * product() const
Definition: ESHandle.h:86
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Definition: value.py:1
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
bool isValid() const
Definition: ESHandle.h:44
unsigned long long uint64_t
Definition: Time.h:13
bool L1TCaloLayer1FetchLUTs(const L1TCaloLayer1FetchLUTsTokens &iTokens, const edm::EventSetup &iSetup, std::vector< std::array< std::array< std::array< uint32_t, nEtBins >, nCalSideBins >, nCalEtaBins > > &eLUT, std::vector< std::array< std::array< std::array< uint32_t, nEtBins >, nCalSideBins >, nCalEtaBins > > &hLUT, std::vector< std::array< std::array< uint32_t, nEtBins >, nHfEtaBins > > &hfLUT, std::vector< unsigned long long int > &hcalFBLUT, std::vector< unsigned int > &ePhiMap, std::vector< unsigned int > &hPhiMap, std::vector< unsigned int > &hfPhiMap, bool useLSB, bool useCalib, bool useECALLUT, bool useHCALLUT, bool useHFLUT, bool useHCALFBLUT, int fwVersion)