CMS 3D CMS Logo

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