CMS 3D CMS Logo

SectorProcessorLUT.cc
Go to the documentation of this file.
2 
3 #include <iostream>
4 #include <fstream>
5 
9 
10 SectorProcessorLUT::SectorProcessorLUT() : version_(0xFFFFFFFF) {}
11 
13 
14 void SectorProcessorLUT::read(bool pc_lut_data, int pc_lut_version) {
15  if (version_ == pc_lut_version)
16  return;
17 
18  edm::LogInfo("L1T") << "EMTF using pc_lut_ver: " << pc_lut_version << ", configured for "
19  << (pc_lut_data ? "data" : "MC");
20 
21  std::string coord_lut_dir = "";
22  if (pc_lut_version == 0)
23  coord_lut_dir = "ph_lut_v1"; // All year 2016
24  else if (pc_lut_version == 1)
25  coord_lut_dir = "ph_lut_v2"; // Beginning of 2017, improved alignment from ideal CMS geometry (MC)
26  else if (pc_lut_version == 2 && pc_lut_data)
27  coord_lut_dir = "ph_lut_v3_data"; // Update in September 2017 from ReReco alignment, data only
28  else if (pc_lut_version == 2)
29  coord_lut_dir = "ph_lut_v2"; // MC still uses ideal CMS aligment
30  else if (pc_lut_version == -1 && pc_lut_data)
31  coord_lut_dir = "ph_lut_v3_data"; // September 2017 data LCT alignment, but use local CPPF LUTs for RPC
32  else if (pc_lut_version == -1)
33  coord_lut_dir = "ph_lut_v2"; // MC using ideal CMS LCT alignment, but use local CPPF LUTs for RPC
34  else
35  throw cms::Exception("L1TMuonEndCap")
36  << "Trying to use EMTF pc_lut_version = " << pc_lut_version << ", does not exist!";
37  // Will catch user trying to run with Global Tag settings on 2016 data, rather than fakeEmtfParams. - AWB 08.06.17
38 
39  std::string coord_lut_path = "L1Trigger/L1TMuon/data/emtf_luts/" + coord_lut_dir + "/";
40 
41  read_file(coord_lut_path + "ph_init_neighbor.txt", ph_init_neighbor_);
42  read_file(coord_lut_path + "ph_disp_neighbor.txt", ph_disp_neighbor_);
43  read_file(coord_lut_path + "th_init_neighbor.txt", th_init_neighbor_);
44  read_file(coord_lut_path + "th_disp_neighbor.txt", th_disp_neighbor_);
45  read_file(coord_lut_path + "th_lut_neighbor.txt", th_lut_neighbor_);
46  read_file(coord_lut_path + "th_corr_lut_neighbor.txt", th_corr_lut_neighbor_);
47 
48  std::string cppf_coord_lut_path = "L1Trigger/L1TMuon/data/cppf/"; // Coordinate LUTs actually used by CPPF
49  bool use_local_cppf_files = (pc_lut_version == -1);
50  if (use_local_cppf_files) { // More accurate coordinate transformation LUTs from Jia Fu
51  cppf_coord_lut_path = "L1Trigger/L1TMuon/data/cppf_luts/angleScale_v1/";
52  }
53 
54  read_cppf_file(cppf_coord_lut_path,
57  use_local_cppf_files); // cppf filenames are hardcoded in the function
58 
59  if (ph_init_neighbor_.size() != 2 * 6 * 61) { // [endcap_2][sector_6][chamber_61]
60  throw cms::Exception("L1TMuonEndCap") << "Expected ph_init_neighbor_ to get " << 2 * 6 * 61 << " values, "
61  << "got " << ph_init_neighbor_.size() << " values.";
62  }
63 
64  if (ph_disp_neighbor_.size() != 2 * 6 * 61) { // [endcap_2][sector_6][chamber_61]
65  throw cms::Exception("L1TMuonEndCap") << "Expected ph_disp_neighbor_ to get " << 2 * 6 * 61 << " values, "
66  << "got " << ph_disp_neighbor_.size() << " values.";
67  }
68 
69  if (th_init_neighbor_.size() != 2 * 6 * 61) { // [endcap_2][sector_6][chamber_61]
70  throw cms::Exception("L1TMuonEndCap") << "Expected th_init_neighbor_ to get " << 2 * 6 * 61 << " values, "
71  << "got " << th_init_neighbor_.size() << " values.";
72  }
73 
74  if (th_disp_neighbor_.size() != 2 * 6 * 61) { // [endcap_2][sector_6][chamber_61]
75  throw cms::Exception("L1TMuonEndCap") << "Expected th_disp_neighbor_ to get " << 2 * 6 * 61 << " values, "
76  << "got " << th_disp_neighbor_.size() << " values.";
77  }
78 
79  if (th_lut_neighbor_.size() != 2 * 6 * 61 * 128) { // [endcap_2][sector_6][chamber_61][wire_128]
80  throw cms::Exception("L1TMuonEndCap") << "Expected th_lut_neighbor_ to get " << 2 * 6 * 61 * 128 << " values, "
81  << "got " << th_lut_neighbor_.size() << " values.";
82  }
83 
84  if (th_corr_lut_neighbor_.size() != 2 * 6 * 7 * 128) { // [endcap_2][sector_6][chamber_61][strip_wire_128]
85  throw cms::Exception("L1TMuonEndCap") << "Expected th_corr_lut_neighbor_ to get " << 2 * 6 * 7 * 128 << " values, "
86  << "got " << th_corr_lut_neighbor_.size() << " values.";
87  }
88 
89  if (cppf_ph_lut_.size() !=
90  2 * 6 * 6 * 6 * 3 *
91  64) { // [endcap_2][rpc_sector_6][rpc_station_ring_6][rpc_subsector_6][rpc_roll_3][rpc_halfstrip_64]
92  throw cms::Exception("L1TMuonEndCap") << "Expected cppf_ph_lut_ to get " << 2 * 6 * 6 * 6 * 3 * 64 << " values, "
93  << "got " << cppf_ph_lut_.size() << " values.";
94  }
95 
96  if (cppf_th_lut_.size() !=
97  2 * 6 * 6 * 6 * 3) { // [endcap_2][rpc_sector_6][rpc_station_ring_6][rpc_subsector_6][rpc_roll_3]
98  throw cms::Exception("L1TMuonEndCap") << "Expected cppf_th_lut_ to get " << 2 * 6 * 6 * 6 * 3 << " values, "
99  << "got " << cppf_th_lut_.size() << " values.";
100  }
101 
102  // clct pattern convertion array from CMSSW
103  //{0.0, 0.0, -0.60, 0.60, -0.64, 0.64, -0.23, 0.23, -0.21, 0.21, 0.0}
104  // 0 0 -5 +5 -5 +5 -2 +2 -2 +2 0
105  ph_patt_corr_ = {0, 0, 5, 5, 5, 5, 2, 2, 2, 2, 0};
106  if (ph_patt_corr_.size() != 11) {
107  throw cms::Exception("L1TMuonEndCap") << "Expected ph_patt_corr_ to get " << 11 << " values, "
108  << "got " << ph_patt_corr_.size() << " values.";
109  }
110 
111  ph_patt_corr_sign_ = {0, 0, 1, 0, 1, 0, 1, 0, 1, 0, 0};
112  if (ph_patt_corr_sign_.size() != 11) {
113  throw cms::Exception("L1TMuonEndCap") << "Expected ph_patt_corr_sign_ to get " << 11 << " values, "
114  << "got " << ph_patt_corr_sign_.size() << " values.";
115  }
116 
117  ph_zone_offset_ = {39, 57, 76, 39, 58, 76, 41, 60, 79, 95, 114, 132, 95, 114, 133, 98, 116, 135,
118  38, 76, 113, 39, 58, 76, 95, 114, 132, 38, 76, 113, 39, 58, 76, 95, 114, 132,
119  38, 76, 113, 38, 57, 76, 95, 113, 132, 21, 21, 23, 1, 21, 1, 21, 1, 20};
120  if (ph_zone_offset_.size() != 6 * 9) {
121  throw cms::Exception("L1TMuonEndCap") << "Expected ph_zone_offset_ to get " << 6 * 9 << " values, "
122  << "got " << ph_zone_offset_.size() << " values.";
123  }
124 
125  // start phi of each chamber in reduced precision, for zone hits,
126  // with negative offset to allow for possible chamber movement
127  ph_init_hard_ = {39, 57, 76, 39, 58, 76, 41, 60, 79, 39, 57, 76, 21, 21, 23, 21, 95, 114, 132, 95,
128  114, 133, 98, 116, 135, 95, 114, 132, 0, 0, 0, 0, 38, 76, 113, 39, 58, 76, 95, 114,
129  132, 1, 21, 0, 0, 0, 0, 0, 38, 76, 113, 39, 58, 76, 95, 114, 132, 1, 21, 0,
130  0, 0, 0, 0, 38, 76, 113, 38, 57, 76, 95, 113, 132, 1, 20, 0, 0, 0, 0, 0};
131  if (ph_init_hard_.size() != 5 * 16) {
132  throw cms::Exception("L1TMuonEndCap") << "Expected ph_init_hard_ to get " << 5 * 16 << " values, "
133  << "got " << ph_init_hard_.size() << " values.";
134  }
135 
136  version_ = pc_lut_version;
137  return;
138 }
139 
140 uint32_t SectorProcessorLUT::get_ph_init(int fw_endcap, int fw_sector, int pc_lut_id) const {
141  const uint32_t index = (fw_endcap * 6 + fw_sector) * 61 + pc_lut_id;
142  uint32_t entry = 0;
143 
144  if (index < ph_init_neighbor_.size()) {
146  } else {
147  edm::LogError("L1T") << "Could not retrieve entry from LUT. fw_endcap: " << fw_endcap
148  << ", fw_sector: " << fw_sector << ", pc_lut_id: " << pc_lut_id;
149  }
150  return entry;
151 }
152 
153 uint32_t SectorProcessorLUT::get_ph_disp(int fw_endcap, int fw_sector, int pc_lut_id) const {
154  const uint32_t index = (fw_endcap * 6 + fw_sector) * 61 + pc_lut_id;
155  uint32_t entry = 0;
156 
157  if (index < ph_disp_neighbor_.size()) {
159  } else {
160  edm::LogError("L1T") << "Could not retrieve entry from LUT. fw_endcap: " << fw_endcap
161  << ", fw_sector: " << fw_sector << ", pc_lut_id: " << pc_lut_id;
162  }
163  return entry;
164 }
165 
166 uint32_t SectorProcessorLUT::get_th_init(int fw_endcap, int fw_sector, int pc_lut_id) const {
167  const uint32_t index = (fw_endcap * 6 + fw_sector) * 61 + pc_lut_id;
168  uint32_t entry = 0;
169 
170  if (index < th_init_neighbor_.size()) {
172  } else {
173  edm::LogError("L1T") << "Could not retrieve entry from LUT. fw_endcap: " << fw_endcap
174  << ", fw_sector: " << fw_sector << ", pc_lut_id: " << pc_lut_id;
175  }
176  return entry;
177 }
178 
179 uint32_t SectorProcessorLUT::get_th_disp(int fw_endcap, int fw_sector, int pc_lut_id) const {
180  const uint32_t index = (fw_endcap * 6 + fw_sector) * 61 + pc_lut_id;
181  uint32_t entry = 0;
182 
183  if (index < th_disp_neighbor_.size()) {
185  } else {
186  edm::LogError("L1T") << "Could not retrieve entry from LUT. fw_endcap: " << fw_endcap
187  << ", fw_sector: " << fw_sector << ", pc_lut_id: " << pc_lut_id;
188  }
189  return entry;
190 }
191 
192 uint32_t SectorProcessorLUT::get_th_lut(int fw_endcap, int fw_sector, int pc_lut_id, int pc_wire_id) const {
193  int pc_lut_id2 = pc_lut_id;
194 
195  // Make ME1/1a the same as ME1/1b
196  if ((9 <= pc_lut_id2 && pc_lut_id2 < 12) || (25 <= pc_lut_id2 && pc_lut_id2 < 28))
197  pc_lut_id2 -= 9;
198  // Make ME1/1a neighbor the same as ME1/1b
199  if (pc_lut_id2 == 15)
200  pc_lut_id2 -= 3;
201 
202  const uint32_t index = ((fw_endcap * 6 + fw_sector) * 61 + pc_lut_id2) * 128 + pc_wire_id;
203  uint32_t entry = 0;
204 
205  if (index < th_lut_neighbor_.size()) {
207  } else {
208  edm::LogError("L1T") << "Could not retrieve entry from LUT. fw_endcap: " << fw_endcap
209  << ", fw_sector: " << fw_sector << ", pc_lut_id: " << pc_lut_id
210  << ", pc_wire_id: " << pc_wire_id;
211  }
212  return entry;
213 }
214 
215 uint32_t SectorProcessorLUT::get_th_corr_lut(int fw_endcap, int fw_sector, int pc_lut_id, int pc_wire_strip_id) const {
216  int pc_lut_id2 = pc_lut_id;
217 
218  // Make ME1/1a the same as ME1/1b
219  if ((9 <= pc_lut_id2 && pc_lut_id2 < 12) || (25 <= pc_lut_id2 && pc_lut_id2 < 28))
220  pc_lut_id2 -= 9;
221  // Make ME1/1a neighbor the same as ME1/1b
222  if (pc_lut_id2 == 15)
223  pc_lut_id2 -= 3;
224 
225  if (pc_lut_id2 <= 3) {
226  pc_lut_id2 -= 0;
227  } else if (pc_lut_id2 == 12) {
228  pc_lut_id2 -= 9;
229  } else if (16 <= pc_lut_id2 && pc_lut_id2 < 19) {
230  pc_lut_id2 -= 12;
231  } else {
232  edm::LogError("L1T") << "get_th_corr_lut(): out of range pc_lut_id: " << pc_lut_id;
233  }
234 
235  const uint32_t index = ((fw_endcap * 6 + fw_sector) * 7 + pc_lut_id2) * 128 + pc_wire_strip_id;
236  uint32_t entry = 0;
237 
238  if (index < th_corr_lut_neighbor_.size()) {
240  } else {
241  edm::LogError("L1T") << "Could not retrieve entry from LUT. fw_endcap: " << fw_endcap
242  << ", fw_sector: " << fw_sector << ", pc_lut_id: " << pc_lut_id
243  << ", pc_wire_strip_id: " << pc_wire_strip_id;
244  }
245  return entry;
246 }
247 
249  const uint32_t index = pattern;
250  uint32_t entry = 0;
251 
252  if (index < ph_patt_corr_.size()) {
253  entry = ph_patt_corr_.at(index);
254  } else {
255  edm::LogError("L1T") << "Could not retrieve entry from LUT. pattern: " << pattern;
256  }
257  return entry;
258 }
259 
261  const uint32_t index = pattern;
262  uint32_t entry = 0;
263 
264  if (index < ph_patt_corr_sign_.size()) {
266  } else {
267  edm::LogError("L1T") << "Could not retrieve entry from LUT. pattern: " << pattern;
268  }
269  return entry;
270 }
271 
272 uint32_t SectorProcessorLUT::get_ph_zone_offset(int pc_station, int pc_chamber) const {
273  const uint32_t index = pc_station * 9 + pc_chamber;
274  uint32_t entry = 0;
275 
276  if (index < ph_zone_offset_.size()) {
278  } else {
279  edm::LogError("L1T") << "Could not retrieve entry from LUT. pc_station: " << pc_station
280  << ", pc_chamber: " << pc_chamber;
281  }
282  return entry;
283 }
284 
285 uint32_t SectorProcessorLUT::get_ph_init_hard(int fw_station, int fw_cscid) const {
286  const uint32_t index = fw_station * 16 + fw_cscid;
287  uint32_t entry = 0;
288 
289  if (index < ph_init_hard_.size()) {
290  entry = ph_init_hard_.at(index);
291  } else {
292  edm::LogError("L1T") << "Could not retrieve entry from LUT. fw_station: " << fw_station
293  << ", fw_cscid: " << fw_cscid;
294  }
295  return entry;
296 }
297 
299  int rpc_region, int rpc_sector, int rpc_station, int rpc_ring, int rpc_subsector, int rpc_roll) const {
300  uint32_t iendcap = (rpc_region == -1) ? 1 : 0;
301  uint32_t isector = (rpc_sector - 1);
302  uint32_t istationring = (rpc_station >= 3) ? ((rpc_station - 3) * 2 + (rpc_ring - 2) + 2) : (rpc_station - 1);
303  uint32_t isubsector = (rpc_subsector - 1);
304  uint32_t iroll = (rpc_roll - 1);
305  return ((((iendcap * 6 + isector) * 6 + istationring) * 6 + isubsector) * 3 + iroll);
306 }
307 
308 uint32_t SectorProcessorLUT::get_cppf_ph_lut(int rpc_region,
309  int rpc_sector,
310  int rpc_station,
311  int rpc_ring,
312  int rpc_subsector,
313  int rpc_roll,
314  int halfstrip,
315  bool is_neighbor) const {
316  const uint32_t th_index = get_cppf_lut_id(rpc_region, rpc_sector, rpc_station, rpc_ring, rpc_subsector, rpc_roll);
317  const uint32_t ph_index = (th_index * 64) + (halfstrip - 1);
318  uint32_t ph = 0;
319 
320  if (ph_index < cppf_ph_lut_.size()) {
321  ph = cppf_ph_lut_.at(ph_index);
322  } else {
323  edm::LogError("L1T") << "Could not retrieve entry from LUT. rpc_region: " << rpc_region
324  << ", rpc_sector: " << rpc_sector << ", rpc_station: " << rpc_station
325  << ", rpc_ring: " << rpc_ring << ", rpc_subsector: " << rpc_subsector
326  << ", rpc_roll: " << rpc_roll << ", halfstrip: " << halfstrip
327  << ", is_neighbor: " << is_neighbor;
328  }
329 
330  if (!is_neighbor && rpc_subsector == 2)
331  ph += 900;
332  return ph;
333 }
334 
336  int rpc_region, int rpc_sector, int rpc_station, int rpc_ring, int rpc_subsector, int rpc_roll) const {
337  const uint32_t th_index = get_cppf_lut_id(rpc_region, rpc_sector, rpc_station, rpc_ring, rpc_subsector, rpc_roll);
338  uint32_t th = 0;
339 
340  if (th_index < cppf_th_lut_.size()) {
341  th = cppf_th_lut_.at(th_index);
342  } else {
343  edm::LogError("L1T") << "Could not retrieve entry from LUT. rpc_region: " << rpc_region
344  << ", rpc_sector: " << rpc_sector << ", rpc_station: " << rpc_station
345  << ", rpc_ring: " << rpc_ring << ", rpc_subsector: " << rpc_subsector
346  << ", rpc_roll: " << rpc_roll;
347  }
348  return th;
349 }
350 
351 void SectorProcessorLUT::read_file(const std::string& filename, std::vector<uint32_t>& vec) {
352  vec.clear();
353 
354  std::ifstream infile;
355  infile.open(edm::FileInPath(filename).fullPath().c_str());
356 
357  int buf;
358  while (infile >> buf) {
359  buf = (buf == -999) ? 0 : buf;
360  vec.push_back(buf);
361  }
362  infile.close();
363 }
364 
366  std::vector<uint32_t>& vec1,
367  std::vector<uint32_t>& vec2,
368  bool local) {
369  auto get_rpc_region = [](uint32_t id) { return (static_cast<int>((id >> 0) & 0X3) + (-1)); };
370  auto get_rpc_sector = [](uint32_t id) { return (static_cast<int>((id >> 7) & 0XF) + (1)); };
371  auto get_rpc_ring = [](uint32_t id) { return (static_cast<int>((id >> 2) & 0X7) + (1)); };
372  auto get_rpc_station = [](uint32_t id) { return (static_cast<int>((id >> 5) & 0X3) + (1)); };
373  auto get_rpc_subsector = [](uint32_t id) { return (static_cast<int>((id >> 12) & 0X7) + (1)); };
374  auto get_rpc_roll = [](uint32_t id) { return (static_cast<int>((id >> 15) & 0X7) + (0)); };
375 
376  std::vector<std::string> cppf_filenames = {
377  "angleScale_RPC_CPPFp1.txt",
378  "angleScale_RPC_CPPFp2.txt",
379  "angleScale_RPC_CPPFp3.txt",
380  "angleScale_RPC_CPPFp4.txt",
381  "angleScale_RPC_CPPFn1.txt",
382  "angleScale_RPC_CPPFn2.txt",
383  "angleScale_RPC_CPPFn3.txt",
384  "angleScale_RPC_CPPFn4.txt",
385  };
386 
387  vec1.clear();
388  vec2.clear();
389  vec1.resize(2 * 6 * 6 * 6 * 3 * 64, 0);
390  vec2.resize(2 * 6 * 6 * 6 * 3, 0);
391 
392  for (size_t i = 0; i < cppf_filenames.size(); ++i) {
393  std::ifstream infile;
394  infile.open(edm::FileInPath(filename + cppf_filenames.at(i)).fullPath().c_str());
395 
396  // std::cout << "\n\nOpening CPPF LUT file " << cppf_filenames.at(i) << std::endl;
397 
398  int buf1, buf2, buf3, buf4, buf5, buf6;
399  // Special variables for transforming centrally-provided CPPF LUTs
400  int buf1_prev = 0, buf2_prev = 0, halfstrip_prev = 0; // Values from previous line in file
401  int line_num = 0; // Line number in file
402  int count_dir = -1; // Direction of half-strip counting: +1 is up, -1 is down
403  int dStrip = 0; // Offset for half-strip from full strip
404  while ((infile >> buf1) && (infile >> buf2) && (infile >> buf3) && (infile >> buf4) && (infile >> buf5) &&
405  (infile >> buf6)) {
406  if ((line_num % 192) == 191)
407  line_num += 1; // Gap in central files vs. Jia Fu's files
408  line_num += 1;
409  // On roughly every-other line, files in L1Trigger/L1TMuon/data/cppf have 0 in the first three columns
410  // Skips a "0 0 0" line once every 192 lines
411  if ((line_num % 2) == 1) {
412  buf1_prev = buf1;
413  buf2_prev = buf2;
414  }
415 
416  if (local && (buf1 == 0 || buf2 == 0)) {
417  throw cms::Exception("L1TMuonEndCap") << "Expected non-0 values, got buf1 = " << buf1 << ", buf2 = " << buf2;
418  }
419  if (!local && (buf1_prev == 0 || buf2_prev == 0)) {
420  throw cms::Exception("L1TMuonEndCap")
421  << "Expected non-0 values, got buf1_prev = " << buf1_prev << ", buf2_prev = " << buf2_prev;
422  }
423 
424  uint32_t id = (local ? buf1 : buf1_prev);
425  int32_t rpc_region = get_rpc_region(id);
426  int32_t rpc_sector = get_rpc_sector(id);
427  int32_t rpc_station = get_rpc_station(id);
428  int32_t rpc_ring = get_rpc_ring(id);
429  int32_t rpc_subsector = get_rpc_subsector(id);
430  int32_t rpc_roll = get_rpc_roll(id);
431 
432  // Offset into halfstrips from centrally-provided LUTs
433  if (buf2_prev * 2 > halfstrip_prev + 8 ||
434  buf2_prev * 2 < halfstrip_prev - 8) { // Starting a new series of strips
435  if (buf2_prev == 1)
436  count_dir = +1; // Starting from a low number, counting up
437  else
438  count_dir = -1; // Starting from a high number, counting down
439  }
440  if (count_dir == -1)
441  dStrip = (buf2_prev * 2 == halfstrip_prev ? 1 : 0);
442  if (count_dir == +1)
443  dStrip = (buf2_prev * 2 == halfstrip_prev + 2 ? 1 : 0);
444  if (buf2_prev * 2 < halfstrip_prev - 8 && buf2_prev == 1)
445  dStrip = 1;
446 
447  //uint32_t strip = buf2;
448  uint32_t halfstrip =
449  (local ? buf2
450  : buf2_prev * 2 -
451  dStrip); // I modified the local text files to use 'halfstrip' instead of 'strip' in column 2
452  halfstrip_prev = halfstrip;
453 
454  uint32_t ph = buf5;
455  uint32_t th = buf6;
456 
457  const uint32_t th_index = get_cppf_lut_id(rpc_region, rpc_sector, rpc_station, rpc_ring, rpc_subsector, rpc_roll);
458  const uint32_t ph_index = (th_index * 64) + (halfstrip - 1);
459 
460  // std::cout << id << " " << rpc_region << " " << rpc_sector << " " << rpc_station << " " << rpc_ring << " "
461  // << rpc_subsector << " " << rpc_roll << " " << halfstrip << " " << th_index << " " << ph_index << std::endl;
462 
463  vec1.at(ph_index) = ph;
464  if (halfstrip == 1)
465  vec2.at(th_index) = th;
466 
467  // Fill gap in centrally-provided LUTs once every 192 lines
468  if (!local && (line_num % 192) == 191)
469  vec1.at(ph_index + 1) = ph;
470 
471  } // End while ((infile >> buf1) && ... && (infile >> buf6))
472  infile.close();
473  } // End loop over CPPF LUT files
474 }
uint32_t get_ph_init_hard(int fw_station, int fw_cscid) const
uint32_t get_ph_patt_corr_sign(int pattern) const
uint32_t get_th_corr_lut(int fw_endcap, int fw_sector, int pc_lut_id, int pc_wire_strip_id) const
uint32_t get_th_lut(int fw_endcap, int fw_sector, int pc_lut_id, int pc_wire_id) const
std::vector< uint32_t > ph_zone_offset_
std::vector< uint32_t > cppf_ph_lut_
void read(bool pc_lut_data, int pc_lut_version)
std::vector< uint32_t > th_corr_lut_neighbor_
std::vector< uint32_t > cppf_th_lut_
Log< level::Error, false > LogError
std::vector< uint32_t > th_disp_neighbor_
uint32_t get_ph_patt_corr(int pattern) const
std::vector< uint32_t > ph_patt_corr_
std::vector< double > vec1
Definition: HCALResponse.h:15
std::vector< uint32_t > ph_init_hard_
uint32_t get_ph_zone_offset(int pc_station, int pc_chamber) const
uint32_t get_cppf_th_lut(int rpc_region, int rpc_sector, int rpc_station, int rpc_ring, int rpc_subsector, int rpc_roll) const
std::vector< uint32_t > ph_patt_corr_sign_
Log< level::Info, false > LogInfo
uint32_t get_cppf_lut_id(int rpc_region, int rpc_sector, int rpc_station, int rpc_ring, int rpc_subsector, int rpc_roll) const
void read_cppf_file(const std::string &filename, std::vector< uint32_t > &vec1, std::vector< uint32_t > &vec2, bool local)
std::vector< uint32_t > th_init_neighbor_
uint32_t get_cppf_ph_lut(int rpc_region, int rpc_sector, int rpc_station, int rpc_ring, int rpc_subsector, int rpc_roll, int halfstrip, bool is_neighbor) const
uint32_t get_ph_disp(int fw_endcap, int fw_sector, int pc_lut_id) const
uint32_t get_ph_init(int fw_endcap, int fw_sector, int pc_lut_id) const
std::vector< uint32_t > ph_disp_neighbor_
uint32_t get_th_init(int fw_endcap, int fw_sector, int pc_lut_id) const
std::vector< uint32_t > ph_init_neighbor_
void read_file(const std::string &filename, std::vector< uint32_t > &vec)
std::vector< vec1 > vec2
Definition: HCALResponse.h:16
uint32_t get_th_disp(int fw_endcap, int fw_sector, int pc_lut_id) const
std::vector< uint32_t > th_lut_neighbor_